HLTrigger
special
plugins
HLTTrackWithHits.h
Go to the documentation of this file.
1
#ifndef HLTrigger_HLTTrackWithHits_H
2
8
// system include files
9
#include <memory>
10
11
// user include files
12
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
13
#include "
FWCore/Framework/interface/EDFilter.h
"
14
15
#include "
FWCore/Framework/interface/Event.h
"
16
#include "
FWCore/Framework/interface/MakerMacros.h
"
17
18
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
19
20
#include "
FWCore/MessageService/interface/MessageLogger.h
"
21
22
#include "
HLTrigger/HLTcore/interface/HLTFilter.h
"
23
#include "
DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h
"
24
#include "
DataFormats/TrackReco/interface/Track.h
"
25
#include "
DataFormats/TrackReco/interface/TrackFwd.h
"
26
27
class
HLTTrackWithHits
:
public
HLTFilter
{
28
public
:
29
explicit
HLTTrackWithHits
(
const
edm::ParameterSet
& iConfig)
30
:
HLTFilter
(iConfig),
31
src_
(iConfig.getParameter<
edm
::
InputTag
>(
"src"
)),
32
minN_
(iConfig.getParameter<
int
>(
"MinN"
)),
33
maxN_
(iConfig.getParameter<
int
>(
"MaxN"
)),
34
MinBPX_
(iConfig.getParameter<
int
>(
"MinBPX"
)),
35
MinFPX_
(iConfig.getParameter<
int
>(
"MinFPX"
)),
36
MinPXL_
(iConfig.getParameter<
int
>(
"MinPXL"
)),
37
MinPT_
(iConfig.getParameter<double>(
"MinPT"
)) {
38
srcToken_
= consumes<reco::TrackCollection>(
src_
);
39
}
40
41
~HLTTrackWithHits
()
override
{}
42
43
static
void
fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions) {
44
edm::ParameterSetDescription
desc;
45
makeHLTFilterDescription
(desc);
46
desc.
add
<
edm::InputTag
>(
"src"
,
edm::InputTag
(
""
));
47
desc.
add
<
int
>(
"MinN"
, 0);
48
desc.
add
<
int
>(
"MaxN"
, 99999);
49
desc.
add
<
int
>(
"MinBPX"
, 0);
50
desc.
add
<
int
>(
"MinFPX"
, 0);
51
desc.
add
<
int
>(
"MinPXL"
, 0);
52
desc.
add
<
double
>(
"MinPT"
, 0.);
53
descriptions.
add
(
"hltTrackWithHits"
, desc);
54
}
55
56
private
:
57
bool
hltFilter
(
edm::Event
&
iEvent
,
58
const
edm::EventSetup
&,
59
trigger::TriggerFilterObjectWithRefs
& filterproduct)
const override
{
60
edm::Handle<reco::TrackCollection>
oHandle;
61
iEvent
.getByToken(
srcToken_
, oHandle);
62
int
s
= oHandle->size();
63
int
count
= 0;
64
for
(
int
i
= 0;
i
!=
s
; ++
i
) {
65
const
reco::Track
&
track
= (*oHandle)[
i
];
66
if
(
track
.pt() <
MinPT_
)
67
continue
;
68
const
reco::HitPattern
&
hits
=
track
.hitPattern();
69
if
(
MinBPX_
> 0 &&
hits
.numberOfValidPixelBarrelHits() >=
MinBPX_
) {
70
++
count
;
71
continue
;
72
}
73
if
(
MinFPX_
> 0 &&
hits
.numberOfValidPixelEndcapHits() >=
MinFPX_
) {
74
++
count
;
75
continue
;
76
}
77
if
(
MinPXL_
> 0 &&
hits
.numberOfValidPixelHits() >=
MinPXL_
) {
78
++
count
;
79
continue
;
80
}
81
}
82
83
bool
answer
= (
count
>=
minN_
&&
count
<=
maxN_
);
84
LogDebug
(
"HLTTrackWithHits"
) <<
module
(
iEvent
) <<
" sees: "
<<
s
<<
" objects. Only: "
<<
count
85
<<
" satisfy the hit requirement. Filter answer is: "
<< (
answer
?
"true"
:
"false"
)
86
<< std::endl;
87
return
answer
;
88
}
89
90
edm::InputTag
src_
;
91
edm::EDGetTokenT<reco::TrackCollection>
srcToken_
;
92
int
minN_
,
maxN_
,
MinBPX_
,
MinFPX_
,
MinPXL_
;
93
double
MinPT_
;
94
};
95
96
#endif
trigger::TriggerFilterObjectWithRefs
Definition:
TriggerFilterObjectWithRefs.h:35
mps_fire.i
i
Definition:
mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition:
ParameterSetDescription.h:95
hfClusterShapes_cfi.hits
hits
Definition:
hfClusterShapes_cfi.py:5
HLTTrackWithHits::HLTTrackWithHits
HLTTrackWithHits(const edm::ParameterSet &iConfig)
Definition:
HLTTrackWithHits.h:29
HLTTrackWithHits::MinFPX_
int MinFPX_
Definition:
HLTTrackWithHits.h:92
edm::EDGetTokenT< reco::TrackCollection >
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::ParameterSetDescription
Definition:
ParameterSetDescription.h:52
EDFilter.h
TriggerFilterObjectWithRefs.h
HLTTrackWithHits::maxN_
int maxN_
Definition:
HLTTrackWithHits.h:92
edm::Handle< reco::TrackCollection >
HLTTrackWithHits
Definition:
HLTTrackWithHits.h:27
HLTFilter
Definition:
HLTFilter.h:28
MessageLogger.h
MakerMacros.h
alignCSCRings.s
s
Definition:
alignCSCRings.py:92
reco::HitPattern
Definition:
HitPattern.h:147
Track.h
HLTTrackWithHits::srcToken_
edm::EDGetTokenT< reco::TrackCollection > srcToken_
Definition:
HLTTrackWithHits.h:91
TrackFwd.h
HLTTrackWithHits::~HLTTrackWithHits
~HLTTrackWithHits() override
Definition:
HLTTrackWithHits.h:41
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition:
ConfigurationDescriptions.cc:57
HLTTrackWithHits::src_
edm::InputTag src_
Definition:
HLTTrackWithHits.h:90
HLTFilter.h
reco::Track
Definition:
Track.h:27
edm::ConfigurationDescriptions
Definition:
ConfigurationDescriptions.h:28
submit.answer
answer
Definition:
submit.py:45
HLTTrackWithHits::MinPT_
double MinPT_
Definition:
HLTTrackWithHits.h:93
HLTTrackWithHits::MinPXL_
int MinPXL_
Definition:
HLTTrackWithHits.h:92
HLT_2018_cff.InputTag
InputTag
Definition:
HLT_2018_cff.py:79016
LogDebug
#define LogDebug(id)
Definition:
MessageLogger.h:670
edm::ParameterSet
Definition:
ParameterSet.h:36
Event.h
KineDebug3::count
void count()
Definition:
KinematicConstrainedVertexUpdatorT.h:21
HLTTrackWithHits::hltFilter
bool hltFilter(edm::Event &iEvent, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition:
HLTTrackWithHits.h:57
HLTTrackWithHits::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition:
HLTTrackWithHits.h:43
createfilelist.int
int
Definition:
createfilelist.py:10
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::EventSetup
Definition:
EventSetup.h:57
HLTTrackWithHits::minN_
int minN_
Definition:
HLTTrackWithHits.h:92
Frameworkfwd.h
HLTFilter::makeHLTFilterDescription
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition:
HLTFilter.cc:25
HLTFilter::module
int module(edm::Event const &) const
Definition:
HLTFilter.cc:47
HLT_2018_cff.track
track
Definition:
HLT_2018_cff.py:10352
ParameterSet.h
edm::Event
Definition:
Event.h:73
HLTTrackWithHits::MinBPX_
int MinBPX_
Definition:
HLTTrackWithHits.h:92
edm::InputTag
Definition:
InputTag.h:15
Generated for CMSSW Reference Manual by
1.8.16