RecoJets
JetAssociationProducers
src
JetSignalVertexCompatibility.cc
Go to the documentation of this file.
1
#include <memory>
2
3
#include "
FWCore/Framework/interface/EDProducer.h
"
4
#include "
FWCore/Framework/interface/Event.h
"
5
#include "
FWCore/Framework/interface/EventSetup.h
"
6
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
7
#include "
FWCore/Utilities/interface/InputTag.h
"
8
9
#include "
DataFormats/Common/interface/Handle.h
"
10
#include "
DataFormats/VertexReco/interface/Vertex.h
"
11
#include "
DataFormats/VertexReco/interface/VertexFwd.h
"
12
#include "
DataFormats/JetReco/interface/Jet.h
"
13
#include "
DataFormats/JetReco/interface/JetTracksAssociation.h
"
14
#include "
DataFormats/JetReco/interface/JetFloatAssociation.h
"
15
16
#include "
TrackingTools/TransientTrack/interface/TransientTrackBuilder.h
"
17
#include "
TrackingTools/Records/interface/TransientTrackRecord.h
"
18
19
#include "
RecoJets/JetAssociationAlgorithms/interface/JetSignalVertexCompatibilityAlgo.h
"
20
21
#include "
JetSignalVertexCompatibility.h
"
22
23
using namespace
reco
;
24
25
JetSignalVertexCompatibility::JetSignalVertexCompatibility
(
const
edm::ParameterSet
&
params
)
26
:
algo
(
params
.getParameter<double>(
"cut"
),
params
.getParameter<double>(
"temperature"
)) {
27
jetTracksAssocToken
= consumes<JetTracksAssociationCollection>(
params
.getParameter<
edm::InputTag
>(
"jetTracksAssoc"
));
28
primaryVerticesToken
= consumes<VertexCollection>(
params
.getParameter<
edm::InputTag
>(
"primaryVertices"
));
29
produces<JetFloatAssociation::Container>();
30
}
31
32
JetSignalVertexCompatibility::~JetSignalVertexCompatibility
() {}
33
34
void
JetSignalVertexCompatibility::produce
(
edm::Event
&
event
,
const
edm::EventSetup
&es) {
35
edm::ESHandle<TransientTrackBuilder>
trackBuilder;
36
es.
get
<
TransientTrackRecord
>().
get
(
"TransientTrackBuilder"
, trackBuilder);
37
38
algo
.
resetEvent
(trackBuilder.
product
());
39
40
edm::Handle<JetTracksAssociationCollection>
jetTracksAssoc
;
41
event
.getByToken(
jetTracksAssocToken
,
jetTracksAssoc
);
42
43
edm::Handle<VertexCollection>
primaryVertices
;
44
event
.getByToken(
primaryVerticesToken
,
primaryVertices
);
45
46
auto
result
= std::make_unique<JetFloatAssociation::Container>(
jetTracksAssoc
->keyProduct());
47
48
for
(
JetTracksAssociationCollection::const_iterator
iter =
jetTracksAssoc
->begin(); iter !=
jetTracksAssoc
->end();
49
++iter) {
50
if
(
primaryVertices
->empty())
51
(*
result
)[iter->first] = -1.;
52
53
const
TrackRefVector
&
tracks
= iter->second;
54
std::vector<float> compatibility =
algo
.
compatibility
(*
primaryVertices
,
tracks
);
55
56
// the first vertex is the presumed signal vertex
57
(*result)[iter->first] = compatibility[0];
58
}
59
60
algo
.
resetEvent
(
nullptr
);
61
62
event
.put(
std::move
(
result
));
63
}
edm::ESHandle::product
T const * product() const
Definition:
ESHandle.h:86
JetTracksAssociation.h
Handle.h
PDWG_EXOHSCP_cff.tracks
tracks
Definition:
PDWG_EXOHSCP_cff.py:28
EDProducer.h
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
JetFloatAssociation.h
reco::JetSignalVertexCompatibilityAlgo::compatibility
std::vector< float > compatibility(const reco::VertexCollection &vertices, const reco::TrackRefVector &tracks) const
Definition:
JetSignalVertexCompatibilityAlgo.cc:63
Jet.h
edm::RefVector< TrackCollection >
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
edm::Handle
Definition:
AssociativeIterator.h:50
JetSignalVertexCompatibility::produce
void produce(edm::Event &event, const edm::EventSetup &es) override
Definition:
JetSignalVertexCompatibility.cc:34
reco::JetSignalVertexCompatibilityAlgo::resetEvent
void resetEvent(const TransientTrackBuilder *trackBuilder)
Definition:
JetSignalVertexCompatibilityAlgo.cc:91
cmsdt::algo
algo
Definition:
constants.h:164
JetSignalVertexCompatibility::primaryVerticesToken
edm::EDGetTokenT< reco::VertexCollection > primaryVerticesToken
Definition:
JetSignalVertexCompatibility.h:25
edm::EventSetup::get
T get() const
Definition:
EventSetup.h:80
TransientTrackRecord
Definition:
TransientTrackRecord.h:11
edm::ESHandle< TransientTrackBuilder >
JetSignalVertexCompatibility.h
Vertex.h
TransientTrackBuilder.h
edm::ParameterSet
Definition:
ParameterSet.h:47
zMuMuMuonUserData.primaryVertices
primaryVertices
Definition:
zMuMuMuonUserData.py:12
Event.h
JetSignalVertexCompatibility::JetSignalVertexCompatibility
JetSignalVertexCompatibility(const edm::ParameterSet ¶ms)
Definition:
JetSignalVertexCompatibility.cc:25
JetSignalVertexCompatibility::~JetSignalVertexCompatibility
~JetSignalVertexCompatibility() override
Definition:
JetSignalVertexCompatibility.cc:32
JetSignalVertexCompatibility::algo
reco::JetSignalVertexCompatibilityAlgo algo
Definition:
JetSignalVertexCompatibility.h:22
edm::AssociationVector::const_iterator
transient_vector_type::const_iterator const_iterator
Definition:
AssociationVector.h:106
edm::EventSetup
Definition:
EventSetup.h:57
TransientTrackRecord.h
get
#define get
InputTag.h
VertexFwd.h
eostools.move
def move(src, dest)
Definition:
eostools.py:511
JetSignalVertexCompatibilityAlgo.h
EventSetup.h
mps_fire.result
result
Definition:
mps_fire.py:311
ParameterSet.h
event
Definition:
event.py:1
edm::Event
Definition:
Event.h:73
JetSignalVertexCompatibility::jetTracksAssocToken
edm::EDGetTokenT< reco::JetTracksAssociationCollection > jetTracksAssocToken
Definition:
JetSignalVertexCompatibility.h:24
edm::InputTag
Definition:
InputTag.h:15
ic5JetVertexCompatibility_cfi.jetTracksAssoc
jetTracksAssoc
Definition:
ic5JetVertexCompatibility_cfi.py:8
Generated for CMSSW Reference Manual by
1.8.16