L1Trigger
L1THGCalUtilities
plugins
selectors
HGC3DClusterSimpleSelector.cc
Go to the documentation of this file.
1
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
2
#include "
FWCore/Framework/interface/global/EDProducer.h
"
3
#include "
FWCore/Framework/interface/Event.h
"
4
#include "
FWCore/Framework/interface/MakerMacros.h
"
5
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
6
7
#include "
DataFormats/L1THGCal/interface/HGCalMulticluster.h
"
8
#include "
CommonTools/Utils/interface/StringCutObjectSelector.h
"
9
10
namespace
l1t
{
11
class
HGC3DClusterSimpleSelector
:
public
edm::global::EDProducer
<> {
12
public
:
13
explicit
HGC3DClusterSimpleSelector
(
const
edm::ParameterSet
&);
14
~HGC3DClusterSimpleSelector
()
override
{}
15
void
produce
(
edm::StreamID
,
edm::Event
&,
const
edm::EventSetup
&)
const override
;
16
17
private
:
18
const
edm::EDGetTokenT<l1t::HGCalMulticlusterBxCollection>
src_
;
19
const
StringCutObjectSelector<l1t::HGCalMulticluster>
cut_
;
20
21
};
// class
22
}
// namespace l1t
23
24
l1t::HGC3DClusterSimpleSelector::HGC3DClusterSimpleSelector
(
const
edm::ParameterSet
&iConfig)
25
: src_(consumes<
l1t
::
HGCalMulticlusterBxCollection
>(iConfig.getParameter<
edm
::
InputTag
>(
"src"
))),
26
cut_(iConfig.getParameter<
std
::
string
>(
"cut"
)) {
27
produces<l1t::HGCalMulticlusterBxCollection>();
28
}
29
30
void
l1t::HGC3DClusterSimpleSelector::produce
(
edm::StreamID
,
edm::Event
&
iEvent
,
const
edm::EventSetup
&)
const
{
31
std::unique_ptr<l1t::HGCalMulticlusterBxCollection>
out
= std::make_unique<l1t::HGCalMulticlusterBxCollection>();
32
33
edm::Handle<l1t::HGCalMulticlusterBxCollection>
multiclusters;
34
iEvent
.getByToken(src_, multiclusters);
35
36
for
(
int
bx
= multiclusters->
getFirstBX
();
bx
<= multiclusters->
getLastBX
(); ++
bx
) {
37
for
(
auto
it = multiclusters->
begin
(
bx
), ed = multiclusters->
end
(
bx
); it != ed; ++it) {
38
const
auto
&
c
= *it;
39
if
(cut_(
c
)) {
40
out
->push_back(
bx
,
c
);
41
}
42
}
43
}
44
45
iEvent
.put(
std::move
(
out
));
46
}
47
using
l1t::HGC3DClusterSimpleSelector
;
48
DEFINE_FWK_MODULE
(
HGC3DClusterSimpleSelector
);
edm::StreamID
Definition:
StreamID.h:30
l1t::HGC3DClusterSimpleSelector::src_
const edm::EDGetTokenT< l1t::HGCalMulticlusterBxCollection > src_
Definition:
HGC3DClusterSimpleSelector.cc:18
edm::EDGetTokenT
Definition:
EDGetToken.h:33
edm
HLT enums.
Definition:
AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition:
HLT_FULL_cff.py:89281
l1GtPatternGenerator_cfi.bx
bx
Definition:
l1GtPatternGenerator_cfi.py:18
l1t::HGC3DClusterSimpleSelector::HGC3DClusterSimpleSelector
HGC3DClusterSimpleSelector(const edm::ParameterSet &)
Definition:
HGC3DClusterSimpleSelector.cc:24
edm::Handle
Definition:
AssociativeIterator.h:50
l1t::HGC3DClusterSimpleSelector::cut_
const StringCutObjectSelector< l1t::HGCalMulticluster > cut_
Definition:
HGC3DClusterSimpleSelector.cc:19
l1t::HGC3DClusterSimpleSelector::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition:
HGC3DClusterSimpleSelector.cc:30
BXVector
Definition:
BXVector.h:15
BXVector::getFirstBX
int getFirstBX() const
MakerMacros.h
HGCalMulticluster.h
l1t::HGC3DClusterSimpleSelector
Definition:
HGC3DClusterSimpleSelector.cc:11
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
l1t::HGC3DClusterSimpleSelector::~HGC3DClusterSimpleSelector
~HGC3DClusterSimpleSelector() override
Definition:
HGC3DClusterSimpleSelector.cc:14
BXVector::begin
const_iterator begin(int bx) const
edm::global::EDProducer
Definition:
EDProducer.h:32
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
BXVector::end
const_iterator end(int bx) const
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
l1t
delete x;
Definition:
CaloConfig.h:22
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::EventSetup
Definition:
EventSetup.h:58
eostools.move
def move(src, dest)
Definition:
eostools.py:511
std
Definition:
JetResolutionObject.h:76
StringCutObjectSelector.h
Frameworkfwd.h
StringCutObjectSelector< l1t::HGCalMulticluster >
MillePedeFileConverter_cfg.out
out
Definition:
MillePedeFileConverter_cfg.py:31
ParameterSet.h
EDProducer.h
c
auto & c
Definition:
CAHitNtupletGeneratorKernelsImpl.h:46
edm::Event
Definition:
Event.h:73
BXVector::getLastBX
int getLastBX() const
Generated for CMSSW Reference Manual by
1.8.16