RecoJets
JetAssociationProducers
src
JetTracksAssociatorAtVertex.cc
Go to the documentation of this file.
1
// \class JetTracksAssociatorAtVertex JetTracksAssociatorAtVertex.cc
2
//
3
// Original Author: Andrea Rizzi
4
// Created: Wed Apr 12 11:12:49 CEST 2006
5
// Accommodated for Jet Package by: Fedor Ratnikov Jul. 30, 2007
6
//
7
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
8
9
// user include files
10
#include "
FWCore/Framework/interface/Event.h
"
11
#include "
DataFormats/Common/interface/View.h
"
12
#include "
DataFormats/JetReco/interface/Jet.h
"
13
#include "
DataFormats/TrackReco/interface/Track.h
"
14
#include "
DataFormats/JetReco/interface/JetTracksAssociation.h
"
15
16
#include "
JetTracksAssociatorAtVertex.h
"
17
18
JetTracksAssociatorAtVertex::JetTracksAssociatorAtVertex
(
const
edm::ParameterSet
& fConfig)
19
: mAssociator(fConfig.getParameter<double>(
"coneSize"
)),
20
mAssociatorAssigned(fConfig.getParameter<double>(
"coneSize"
)),
21
useAssigned
(
false
),
22
pvSrc
() {
23
mJets
= consumes<edm::View<reco::Jet> >(fConfig.
getParameter
<
edm::InputTag
>(
"jets"
));
24
mTracks
= consumes<reco::TrackCollection>(fConfig.
getParameter
<
edm::InputTag
>(
"tracks"
));
25
if
(fConfig.
exists
(
"useAssigned"
)) {
26
useAssigned
= fConfig.
getParameter
<
bool
>(
"useAssigned"
);
27
pvSrc
= consumes<reco::VertexCollection>(fConfig.
getParameter
<
edm::InputTag
>(
"pvSrc"
));
28
}
29
30
produces<reco::JetTracksAssociation::Container>();
31
}
32
33
JetTracksAssociatorAtVertex::~JetTracksAssociatorAtVertex
() {}
34
35
void
JetTracksAssociatorAtVertex::produce
(
edm::Event
&
fEvent
,
const
edm::EventSetup
& fSetup) {
36
edm::Handle<edm::View<reco::Jet>
> jets_h;
37
fEvent
.getByToken(
mJets
, jets_h);
38
edm::Handle<reco::TrackCollection>
tracks_h;
39
fEvent
.getByToken(
mTracks
, tracks_h);
40
41
auto
jetTracks
= std::make_unique<reco::JetTracksAssociation::Container>(
reco::JetRefBaseProd
(jets_h));
42
43
// format inputs
44
std::vector<edm::RefToBase<reco::Jet> >
allJets
;
45
allJets
.reserve(jets_h->size());
46
for
(
unsigned
i
= 0;
i
< jets_h->size(); ++
i
)
47
allJets
.push_back(jets_h->refAt(
i
));
48
std::vector<reco::TrackRef>
allTracks
;
49
allTracks
.reserve(tracks_h->size());
50
// run algo
51
for
(
unsigned
i
= 0;
i
< tracks_h->size(); ++
i
) {
52
allTracks
.push_back(
reco::TrackRef
(tracks_h,
i
));
53
}
54
if
(!
useAssigned
) {
55
mAssociator
.
produce
(&*
jetTracks
,
allJets
,
allTracks
);
56
}
else
{
57
edm::Handle<reco::VertexCollection>
pvHandle;
58
fEvent
.getByToken(
pvSrc
, pvHandle);
59
const
reco::VertexCollection
&
vertices
= *pvHandle.
product
();
60
mAssociatorAssigned
.
produce
(&*
jetTracks
,
allJets
,
allTracks
,
vertices
);
61
}
62
63
// store output
64
fEvent
.put(
std::move
(
jetTracks
));
65
}
JetTracksAssociation.h
mps_fire.i
i
Definition:
mps_fire.py:428
MessageLogger.h
funct::false
false
Definition:
Factorize.h:29
JetTracksAssociatorAtVertex::mJets
edm::EDGetTokenT< edm::View< reco::Jet > > mJets
Definition:
JetTracksAssociatorAtVertex.h:27
edm::Handle::product
T const * product() const
Definition:
Handle.h:70
JetTracksAssociatorAtVertex.h
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition:
VertexFwd.h:9
muonTagProbeFilters_cff.allTracks
allTracks
Definition:
muonTagProbeFilters_cff.py:22
Jet.h
edm::Handle
Definition:
AssociativeIterator.h:50
edm::Ref< TrackCollection >
JetTracksAssociatorAtVertex::JetTracksAssociatorAtVertex
JetTracksAssociatorAtVertex(const edm::ParameterSet &)
Definition:
JetTracksAssociatorAtVertex.cc:18
Track.h
JetTracksAssociatorAtVertex::mAssociator
JetTracksAssociationDRVertex mAssociator
Definition:
JetTracksAssociatorAtVertex.h:31
JetTracksAssociatorAtVertex::mAssociatorAssigned
JetTracksAssociationDRVertexAssigned mAssociatorAssigned
Definition:
JetTracksAssociatorAtVertex.h:32
edm::ParameterSet::exists
bool exists(std::string const ¶meterName) const
checks if a parameter exists
Definition:
ParameterSet.cc:681
HLT_FULL_cff.useAssigned
useAssigned
Definition:
HLT_FULL_cff.py:50852
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
HLT_FULL_cff.jetTracks
jetTracks
Definition:
HLT_FULL_cff.py:50866
JetTracksAssociationDRVertex::produce
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks) const
Definition:
JetTracksAssociationDRVertex.cc:14
JetTracksAssociatorAtVertex::~JetTracksAssociatorAtVertex
~JetTracksAssociatorAtVertex() override
Definition:
JetTracksAssociatorAtVertex.cc:33
reco::JetExtendedAssociation::allJets
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
Definition:
JetExtendedAssociation.cc:60
edm::EventSetup
Definition:
EventSetup.h:58
JetTracksAssociatorAtVertex::pvSrc
edm::EDGetTokenT< reco::VertexCollection > pvSrc
if true, use the track/jet association with vertex assignment to tracks
Definition:
JetTracksAssociatorAtVertex.h:34
hcaldqm::fEvent
Definition:
DQTask.h:32
JetTracksAssociatorAtVertex::mTracks
edm::EDGetTokenT< reco::TrackCollection > mTracks
Definition:
JetTracksAssociatorAtVertex.h:28
eostools.move
def move(src, dest)
Definition:
eostools.py:511
HLT_FULL_cff.pvSrc
pvSrc
Definition:
HLT_FULL_cff.py:34856
JetTracksAssociatorAtVertex::useAssigned
bool useAssigned
Definition:
JetTracksAssociatorAtVertex.h:33
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
View.h
edm::Event
Definition:
Event.h:73
JetTracksAssociationDRVertexAssigned::produce
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks, const reco::VertexCollection &vertices) const
Definition:
JetTracksAssociationDRVertexAssigned.cc:16
edm::InputTag
Definition:
InputTag.h:15
JetTracksAssociatorAtVertex::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition:
JetTracksAssociatorAtVertex.cc:35
edm::RefToBaseProd
Definition:
RefToBase.h:65
pwdgSkimBPark_cfi.vertices
vertices
Definition:
pwdgSkimBPark_cfi.py:7
Generated for CMSSW Reference Manual by
1.8.16