src
JetMETCorrections
JetCorrector
plugins
JetTracksAssociationToTrackRefs.cc
Go to the documentation of this file.
1
#include "
FWCore/Framework/interface/global/EDProducer.h
"
2
#include "
FWCore/Framework/interface/MakerMacros.h
"
3
#include "
FWCore/Framework/interface/Event.h
"
4
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
5
#include "
FWCore/ParameterSet/interface/ConfigurationDescriptions.h
"
6
#include "
FWCore/ParameterSet/interface/ParameterSetDescription.h
"
7
8
#include "
DataFormats/Common/interface/EDProductfwd.h
"
9
#include "
DataFormats/TrackReco/interface/Track.h
"
10
#include "
DataFormats/TrackReco/interface/TrackFwd.h
"
11
#include "
DataFormats/JetReco/interface/Jet.h
"
12
#include "
DataFormats/JetReco/interface/JetTracksAssociation.h
"
13
14
#include "
JetMETCorrections/JetCorrector/interface/JetCorrector.h
"
15
16
#include <unordered_set>
17
23
class
JetTracksAssociationToTrackRefs
:
public
edm::global::EDProducer
<> {
24
public
:
25
JetTracksAssociationToTrackRefs
(
const
edm::ParameterSet
& iConfig);
26
27
static
void
fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions);
28
29
void
produce
(
edm::StreamID
,
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup)
const override
;
30
31
private
:
32
edm::EDGetTokenT<reco::JetTracksAssociation::Container>
associationToken_
;
33
edm::EDGetTokenT<edm::View<reco::Jet>
>
jetToken_
;
34
edm::EDGetTokenT<reco::JetCorrector>
correctorToken_
;
35
const
double
ptMin_
;
36
};
37
38
JetTracksAssociationToTrackRefs::JetTracksAssociationToTrackRefs
(
const
edm::ParameterSet
& iConfig)
39
: associationToken_(
40
consumes<
reco
::
JetTracksAssociation
::
Container
>(iConfig.getParameter<
edm
::
InputTag
>(
"association"
))),
41
jetToken_(consumes<
edm
::
View
<
reco
::
Jet
>>(iConfig.getParameter<
edm
::
InputTag
>(
"jets"
))),
42
correctorToken_(consumes<
reco
::JetCorrector>(iConfig.getParameter<
edm
::
InputTag
>(
"corrector"
))),
43
ptMin_(iConfig.getParameter<double>(
"correctedPtMin"
)) {
44
produces<reco::TrackRefVector>();
45
}
46
47
void
JetTracksAssociationToTrackRefs::fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions) {
48
edm::ParameterSetDescription
desc
;
49
desc
.add<
edm::InputTag
>(
"association"
,
edm::InputTag
(
"ak4JetTracksAssociatorAtVertexPF"
));
50
desc
.add<
edm::InputTag
>(
"jets"
,
edm::InputTag
(
"ak4PFJetsCHS"
));
51
desc
.add<
edm::InputTag
>(
"corrector"
,
edm::InputTag
(
"ak4PFL1FastL2L3Corrector"
));
52
desc
.add<
double
>(
"correctedPtMin"
, 0);
53
descriptions.
add
(
"jetTracksAssociationToTrackRefs"
,
desc
);
54
}
55
56
void
JetTracksAssociationToTrackRefs::produce
(
edm::StreamID
,
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup)
const
{
57
edm::Handle<reco::JetTracksAssociation::Container>
h_assoc;
58
iEvent
.getByToken(
associationToken_
, h_assoc);
59
const
reco::JetTracksAssociation::Container
&
association
= *h_assoc;
60
61
edm::Handle<edm::View<reco::Jet>
> h_jets;
62
iEvent
.getByToken(
jetToken_
, h_jets);
63
const
edm::View<reco::Jet>
&
jets
= *h_jets;
64
65
edm::Handle<reco::JetCorrector>
h_corrector;
66
iEvent
.getByToken(
correctorToken_
, h_corrector);
67
const
reco::JetCorrector
&
corrector
= *h_corrector;
68
69
auto
ret
= std::make_unique<reco::TrackRefVector>();
70
std::unordered_set<reco::TrackRefVector::key_type> alreadyAdded;
71
72
// Exctract tracks only for jets passing certain pT threshold after
73
// correction
74
for
(
size_t
i
= 0;
i
<
jets
.size(); ++
i
) {
75
edm::RefToBase<reco::Jet>
jetRef =
jets
.refAt(
i
);
76
const
reco::Jet
&
jet
= *jetRef;
77
78
auto
p4 =
jet
.p4();
79
80
// Energy correction in the most general way
81
if
(!
corrector
.vectorialCorrection()) {
82
double
scale
= 1;
83
if
(!
corrector
.refRequired()) {
84
scale
=
corrector
.correction(
jet
);
85
}
else
{
86
scale
=
corrector
.correction(
jet
, jetRef);
87
}
88
p4 = p4 *
scale
;
89
}
else
{
90
corrector
.correction(
jet
, jetRef, p4);
91
}
92
93
if
(p4.pt() <=
ptMin_
)
94
continue
;
95
96
for
(
const
auto
& trackRef :
association
[jetRef]) {
97
if
(alreadyAdded.find(trackRef.key()) == alreadyAdded.end()) {
98
ret
->push_back(trackRef);
99
alreadyAdded.insert(trackRef.key());
100
}
101
}
102
}
103
104
iEvent
.put(
std::move
(
ret
));
105
}
106
107
DEFINE_FWK_MODULE
(
JetTracksAssociationToTrackRefs
);
metsig::jet
Definition:
SignAlgoResolutions.h:47
JetCorrector.h
JetTracksAssociationToTrackRefs::jetToken_
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition:
JetTracksAssociationToTrackRefs.cc:33
JetTracksAssociation.h
ProducerED_cfi.InputTag
InputTag
Definition:
ProducerED_cfi.py:5
JetTracksAssociationToTrackRefs::produce
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
Definition:
JetTracksAssociationToTrackRefs.cc:56
mps_fire.i
i
Definition:
mps_fire.py:429
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition:
runTheMatrix.py:732
reco::Jet
Base class for all types of Jets.
Definition:
Jet.h:20
edm::StreamID
Definition:
StreamID.h:30
Event.h
JetTracksAssociationToTrackRefs::ptMin_
const double ptMin_
Definition:
JetTracksAssociationToTrackRefs.cc:35
edm::Handle
Definition:
AssociativeIterator.h:50
PDWG_EXODelayedJetMET_cff.jets
jets
Definition:
PDWG_EXODelayedJetMET_cff.py:14
TrackFwd.h
JetTracksAssociationToTrackRefs::JetTracksAssociationToTrackRefs
JetTracksAssociationToTrackRefs(const edm::ParameterSet &iConfig)
Definition:
JetTracksAssociationToTrackRefs.cc:38
JetTracksAssociationToTrackRefs::correctorToken_
edm::EDGetTokenT< reco::JetCorrector > correctorToken_
Definition:
JetTracksAssociationToTrackRefs.cc:34
edm::EDGetTokenT
Definition:
EDGetToken.h:37
edm::RefToBase< reco::Jet >
edm::ParameterSetDescription
Definition:
ParameterSetDescription.h:52
ParameterSet.h
edm::View
Definition:
CaloClusterFwd.h:14
iEvent
int iEvent
Definition:
GenABIO.cc:224
ParameterSetDescription.h
hgcal::association
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
Definition:
LCToCPAssociatorByEnergyScoreImpl.h:62
Jet
Definition:
Jet.py:1
submitPVResolutionJobs.desc
string desc
Definition:
submitPVResolutionJobs.py:251
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
edm::EventSetup
Definition:
EventSetup.h:56
Jet.h
edm::global::EDProducer
Definition:
EDProducer.h:32
pfClustersFromHGC3DClusters_cfi.corrector
corrector
Definition:
pfClustersFromHGC3DClusters_cfi.py:4
reco::JetCorrector
Definition:
JetCorrector.h:33
JetTracksAssociationToTrackRefs::associationToken_
edm::EDGetTokenT< reco::JetTracksAssociation::Container > associationToken_
Definition:
JetTracksAssociationToTrackRefs.cc:32
edm::AssociationVector
Definition:
AssociationVector.h:67
EDProductfwd.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition:
ConfigurationDescriptions.cc:57
JetTracksAssociationToTrackRefs
Definition:
JetTracksAssociationToTrackRefs.cc:23
EDProducer.h
JetTracksAssociation
Association between jets and float value.
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::InputTag
Definition:
InputTag.h:15
pfClustersFromCombinedCaloHF_cfi.scale
scale
Definition:
pfClustersFromCombinedCaloHF_cfi.py:51
Track.h
edm::ParameterSet
Definition:
ParameterSet.h:48
ConfigurationDescriptions.h
JetTracksAssociationToTrackRefs::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition:
JetTracksAssociationToTrackRefs.cc:47
edm::Event
Definition:
Event.h:73
MakerMacros.h
reco::JetExtendedAssociation::Container
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
Definition:
JetExtendedAssociation.h:29
sistrip::View
View
Definition:
ConstantsForView.h:26
eostools.move
def move(src, dest)
Definition:
eostools.py:511
edm::ConfigurationDescriptions
Definition:
ConfigurationDescriptions.h:28
Generated for CMSSW Reference Manual by
1.8.14