Validation
RecoTrack
plugins
CrossingFramePSimHitToPSimHits.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 "
SimDataFormats/TrackingHit/interface/PSimHit.h
"
9
#include "
SimDataFormats/CrossingFrame/interface/CrossingFrame.h
"
10
11
class
CrossingFramePSimHitToPSimHitsConverter
:
public
edm::global::EDProducer
<> {
12
public
:
13
CrossingFramePSimHitToPSimHitsConverter
(
const
edm::ParameterSet
& iConfig);
14
15
static
void
fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions);
16
17
void
produce
(
edm::StreamID
,
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup)
const override
;
18
19
private
:
20
struct
InputInfo
{
21
using
Token
=
edm::EDGetTokenT<CrossingFrame<PSimHit>
>;
22
InputInfo
(
Token
t
,
const
std::string
&
i
) :
token
(
t
),
instance
(
i
) {}
23
Token
token
;
24
std::string
instance
;
25
};
26
27
std::vector<InputInfo>
input_
;
28
};
29
30
CrossingFramePSimHitToPSimHitsConverter::CrossingFramePSimHitToPSimHitsConverter
(
const
edm::ParameterSet
& iConfig) {
31
auto
src
= iConfig.
getParameter
<std::vector<edm::InputTag>>(
"src"
);
32
input_
.reserve(
src
.size());
33
for
(
const
auto
&
tag
:
src
) {
34
input_
.emplace_back(
consumes
<
CrossingFrame<PSimHit>
>(
tag
),
tag
.instance());
35
produces<std::vector<PSimHit>>(
input_
.back().instance);
36
}
37
}
38
39
void
CrossingFramePSimHitToPSimHitsConverter::fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions) {
40
edm::ParameterSetDescription
desc
;
41
desc
.add<std::vector<edm::InputTag>>(
"src"
, std::vector<edm::InputTag>());
42
descriptions.
add
(
"crossingFramePSimHitToPSimHits"
,
desc
);
43
}
44
45
void
CrossingFramePSimHitToPSimHitsConverter::produce
(
edm::StreamID
,
46
edm::Event
&
iEvent
,
47
const
edm::EventSetup
& iSetup)
const
{
48
for
(
const
auto
&
input
:
input_
) {
49
edm::Handle<CrossingFrame<PSimHit>
> hframe;
50
iEvent
.getByToken(
input
.token, hframe);
51
const
auto
&
frame
= *hframe;
52
const
auto
& signalHits =
frame
.getSignal();
53
const
auto
& pileupHits =
frame
.getPileups();
54
55
auto
output
= std::make_unique<std::vector<PSimHit>>();
56
output
->reserve(signalHits.size() + pileupHits.size());
57
for
(
const
auto
& ptr : signalHits)
58
output
->emplace_back(*ptr);
59
for
(
const
auto
& ptr : pileupHits)
60
output
->emplace_back(*ptr);
61
iEvent
.put(
std::move
(
output
),
input
.instance);
62
}
63
}
64
65
//define this as a plug-in
66
DEFINE_FWK_MODULE
(
CrossingFramePSimHitToPSimHitsConverter
);
ConfigurationDescriptions.h
edm::StreamID
Definition:
StreamID.h:30
mps_fire.i
i
Definition:
mps_fire.py:428
input
static const std::string input
Definition:
EdmProvDump.cc:48
CrossingFramePSimHitToPSimHitsConverter::InputInfo::instance
std::string instance
Definition:
CrossingFramePSimHitToPSimHits.cc:24
convertSQLitetoXML_cfg.output
output
Definition:
convertSQLitetoXML_cfg.py:72
edm::EDGetTokenT< CrossingFrame< PSimHit > >
CrossingFramePSimHitToPSimHitsConverter::InputInfo::token
Token token
Definition:
CrossingFramePSimHitToPSimHits.cc:23
CrossingFrame.h
CrossingFramePSimHitToPSimHitsConverter
Definition:
CrossingFramePSimHitToPSimHits.cc:11
edm::ParameterSetDescription
Definition:
ParameterSetDescription.h:52
CrossingFramePSimHitToPSimHitsConverter::input_
std::vector< InputInfo > input_
Definition:
CrossingFramePSimHitToPSimHits.cc:27
edm::Handle
Definition:
AssociativeIterator.h:50
CrossingFrame
Definition:
CrossingFrame.h:37
MakerMacros.h
PSimHit.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition:
ConfigurationDescriptions.cc:57
CrossingFramePSimHitToPSimHitsConverter::produce
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
Definition:
CrossingFramePSimHitToPSimHits.cc:45
CrossingFramePSimHitToPSimHitsConverter::CrossingFramePSimHitToPSimHitsConverter
CrossingFramePSimHitToPSimHitsConverter(const edm::ParameterSet &iConfig)
Definition:
CrossingFramePSimHitToPSimHits.cc:30
ParameterSetDescription.h
CrossingFramePSimHitToPSimHitsConverter::InputInfo::InputInfo
InputInfo(Token t, const std::string &i)
Definition:
CrossingFramePSimHitToPSimHits.cc:22
edm::global::EDProducer
Definition:
EDProducer.h:32
edm::ConfigurationDescriptions
Definition:
ConfigurationDescriptions.h:28
edm::ParameterSet
Definition:
ParameterSet.h:47
TrackRefitter_38T_cff.src
src
Definition:
TrackRefitter_38T_cff.py:24
Event.h
makeGlobalPositionRcd_cfg.tag
tag
Definition:
makeGlobalPositionRcd_cfg.py:6
CrossingFramePSimHitToPSimHitsConverter::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition:
CrossingFramePSimHitToPSimHits.cc:39
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::EventSetup
Definition:
EventSetup.h:58
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
submitPVResolutionJobs.desc
string desc
Definition:
submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition:
eostools.py:511
amptDefault_cfi.frame
frame
Definition:
amptDefault_cfi.py:12
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
ParameterSet.h
edm::EDConsumerBase::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition:
EDConsumerBase.h:159
EDProducer.h
edm::Event
Definition:
Event.h:73
submitPVValidationJobs.t
string t
Definition:
submitPVValidationJobs.py:644
CrossingFramePSimHitToPSimHitsConverter::InputInfo
Definition:
CrossingFramePSimHitToPSimHits.cc:20
Generated for CMSSW Reference Manual by
1.8.16