RecoParticleFlow
PFProducer
plugins
importers
TrackTimingImporter.cc
Go to the documentation of this file.
1
#include "
RecoParticleFlow/PFProducer/interface/BlockElementImporterBase.h
"
2
#include "
DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h
"
3
#include "
DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h
"
4
#include "
DataFormats/ParticleFlowReco/interface/PFRecTrack.h
"
5
#include "
DataFormats/TrackReco/interface/Track.h
"
6
#include "
DataFormats/MuonReco/interface/Muon.h
"
7
#include "
DataFormats/Common/interface/ValueMap.h
"
8
#include "
RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h
"
9
#include "
RecoParticleFlow/PFTracking/interface/PFTrackAlgoTools.h
"
10
11
// this doesn't actually import anything,
12
// but rather applies time stamps to tracks after they are all inserted
13
14
class
TrackTimingImporter
:
public
BlockElementImporterBase
{
15
public
:
16
TrackTimingImporter
(
const
edm::ParameterSet
& conf,
edm::ConsumesCollector
& sumes)
17
:
BlockElementImporterBase
(conf, sumes),
18
useTimeQuality_
(conf.existsAs<
edm
::
InputTag
>(
"timeQualityMap"
)),
19
timeQualityThreshold_
(
useTimeQuality_
? conf.getParameter<double>(
"timeQualityThreshold"
) : -99.),
20
srcTime_
(sumes.consumes<
edm
::ValueMap<
float
>>(conf.getParameter<
edm
::
InputTag
>(
"timeValueMap"
))),
21
srcTimeError_
(sumes.consumes<
edm
::ValueMap<
float
>>(conf.getParameter<
edm
::
InputTag
>(
"timeErrorMap"
))),
22
srcTimeQuality_
(
useTimeQuality_
23
? sumes.consumes<
edm
::ValueMap<
float
>>(conf.getParameter<
edm
::
InputTag
>(
"timeQualityMap"
))
24
:
edm
::EDGetTokenT<
edm
::ValueMap<
float
>>()),
25
srcTimeGsf_
(sumes.consumes<
edm
::ValueMap<
float
>>(conf.getParameter<
edm
::
InputTag
>(
"timeValueMapGsf"
))),
26
srcTimeErrorGsf_
(sumes.consumes<
edm
::ValueMap<
float
>>(conf.getParameter<
edm
::
InputTag
>(
"timeErrorMapGsf"
))),
27
srcTimeQualityGsf_
(
useTimeQuality_
? sumes.consumes<
edm
::ValueMap<
float
>>(
28
conf.getParameter<
edm
::
InputTag
>(
"timeQualityMapGsf"
))
29
:
edm
::EDGetTokenT<
edm
::ValueMap<
float
>>()),
30
debug_
(conf.getUntrackedParameter<
bool
>(
"debug"
,
false
)) {}
31
32
void
importToBlock
(
const
edm::Event
&,
ElementList
&)
const override
;
33
34
private
:
35
const
bool
useTimeQuality_
;
36
const
double
timeQualityThreshold_
;
37
edm::EDGetTokenT<edm::ValueMap<float>
>
srcTime_
,
srcTimeError_
,
srcTimeQuality_
,
srcTimeGsf_
,
srcTimeErrorGsf_
,
38
srcTimeQualityGsf_
;
39
const
bool
debug_
;
40
};
41
42
DEFINE_EDM_PLUGIN
(
BlockElementImporterFactory
,
TrackTimingImporter
,
"TrackTimingImporter"
);
43
44
void
TrackTimingImporter::importToBlock
(
const
edm::Event
&
e
,
BlockElementImporterBase::ElementList
& elems)
const
{
45
typedef
BlockElementImporterBase::ElementList::value_type
ElementType;
46
47
auto
const
&
time
=
e
.get(
srcTime_
);
48
auto
const
& timeErr =
e
.get(
srcTimeError_
);
49
auto
const
& timeGsf =
e
.get(
srcTimeGsf_
);
50
auto
const
& timeErrGsf =
e
.get(
srcTimeErrorGsf_
);
51
52
edm::Handle<edm::ValueMap<float>
> timeQualH, timeQualGsfH;
53
if
(
useTimeQuality_
) {
54
e
.getByToken(
srcTimeQuality_
, timeQualH);
55
e
.getByToken(
srcTimeQualityGsf_
, timeQualGsfH);
56
}
57
58
for
(
auto
& elem : elems) {
59
if
(
reco::PFBlockElement::TRACK
== elem->type()) {
60
const
auto
& ref = elem->trackRef();
61
if
(
time
.contains(ref.id())) {
62
const
bool
assocQuality =
useTimeQuality_
? (*timeQualH)[ref] >
timeQualityThreshold_
:
true
;
63
if
(assocQuality) {
64
elem->setTime(
time
[ref], timeErr[ref]);
65
}
else
{
66
elem->setTime(0., -1.);
67
}
68
}
69
if
(
debug_
) {
70
edm::LogInfo
(
"TrackTimingImporter"
)
71
<<
"Track with pT / eta "
<< ref->pt() <<
" / "
<< ref->eta() <<
" has time: "
<< elem->time() <<
" +/- "
72
<< elem->timeError() << std::endl;
73
}
74
}
else
if
(
reco::PFBlockElement::GSF
== elem->type()) {
75
const
auto
& ref = static_cast<const reco::PFBlockElementGsfTrack*>(elem.get())->GsftrackRef();
76
if
(timeGsf.contains(ref.id())) {
77
const
bool
assocQuality =
useTimeQuality_
? (*timeQualGsfH)[ref] >
timeQualityThreshold_
:
true
;
78
if
(assocQuality) {
79
elem->setTime(timeGsf[ref], timeErrGsf[ref]);
80
}
else
{
81
elem->setTime(0., -1.);
82
}
83
}
84
if
(
debug_
) {
85
edm::LogInfo
(
"TrackTimingImporter"
)
86
<<
"Track with pT / eta "
<< ref->pt() <<
" / "
<< ref->eta() <<
" has time: "
<< elem->time() <<
" +/- "
87
<< elem->timeError() << std::endl;
88
}
89
}
90
}
91
}
electrons_cff.bool
bool
Definition:
electrons_cff.py:393
Muon.h
dqmMemoryStats.float
float
Definition:
dqmMemoryStats.py:127
PFMuonAlgo.h
funct::false
false
Definition:
Factorize.h:29
edm::EDGetTokenT
Definition:
EDGetToken.h:33
edm
HLT enums.
Definition:
AlignableModifier.h:19
TrackTimingImporter::srcTimeError_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeError_
Definition:
TrackTimingImporter.cc:37
HLT_FULL_cff.InputTag
InputTag
Definition:
HLT_FULL_cff.py:85964
TrackTimingImporter::srcTimeGsf_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeGsf_
Definition:
TrackTimingImporter.cc:37
TrackTimingImporter::importToBlock
void importToBlock(const edm::Event &, ElementList &) const override
Definition:
TrackTimingImporter.cc:44
edm::LogInfo
Log< level::Info, false > LogInfo
Definition:
MessageLogger.h:125
edm::Handle
Definition:
AssociativeIterator.h:50
TrackTimingImporter::debug_
const bool debug_
Definition:
TrackTimingImporter.cc:39
Track.h
BlockElementImporterBase.h
reco::PFBlockElement::TRACK
Definition:
PFBlockElement.h:32
TrackTimingImporter
Definition:
TrackTimingImporter.cc:14
TrackTimingImporter::srcTimeErrorGsf_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeErrorGsf_
Definition:
TrackTimingImporter.cc:37
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
BlockElementImporterBase::ElementList
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
Definition:
BlockElementImporterBase.h:16
PFTrackAlgoTools.h
edm::ParameterSet
Definition:
ParameterSet.h:47
reco::PFBlockElement::GSF
Definition:
PFBlockElement.h:37
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
TrackTimingImporter::srcTimeQuality_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeQuality_
Definition:
TrackTimingImporter.cc:37
TrackTimingImporter::srcTime_
edm::EDGetTokenT< edm::ValueMap< float > > srcTime_
Definition:
TrackTimingImporter.cc:37
TrackTimingImporter::timeQualityThreshold_
const double timeQualityThreshold_
Definition:
TrackTimingImporter.cc:36
reco::JetExtendedAssociation::value_type
Container::value_type value_type
Definition:
JetExtendedAssociation.h:30
PFRecTrack.h
ValueMap.h
ntuplemaker.time
time
Definition:
ntuplemaker.py:310
edm::Event
Definition:
Event.h:73
TrackTimingImporter::TrackTimingImporter
TrackTimingImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Definition:
TrackTimingImporter.cc:16
PFBlockElementGsfTrack.h
TrackTimingImporter::srcTimeQualityGsf_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeQualityGsf_
Definition:
TrackTimingImporter.cc:37
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
PFBlockElementTrack.h
MillePedeFileConverter_cfg.e
e
Definition:
MillePedeFileConverter_cfg.py:37
BlockElementImporterBase
Definition:
BlockElementImporterBase.h:14
TrackTimingImporter::useTimeQuality_
const bool useTimeQuality_
Definition:
TrackTimingImporter.cc:35
Generated for CMSSW Reference Manual by
1.8.16