PhysicsTools
FWLite
interface
WSelectorFast.h
Go to the documentation of this file.
1
#ifndef PhysicsTools_FWLite_WSelectorFast_h
2
#define PhysicsTools_FWLite_WSelectorFast_h
3
#include "
DataFormats/PatCandidates/interface/MET.h
"
4
#include "
DataFormats/PatCandidates/interface/Muon.h
"
5
#include "
PhysicsTools/SelectorUtils/interface/EventSelector.h
"
6
17
class
WSelector
:
public
EventSelector
{
18
19
public
:
21
WSelector
(
edm::ParameterSet
const
&
params
) :
22
muonSrc_
(
params
.getParameter<
edm
::
InputTag
>(
"muonSrc"
)),
23
metSrc_
(
params
.getParameter<
edm
::
InputTag
>(
"metSrc"
))
24
{
25
double
muonPtMin
=
params
.getParameter<
double
>(
"muonPtMin"
);
26
double
metMin
=
params
.getParameter<
double
>(
"metMin"
);
27
push_back(
"Muon Pt"
,
muonPtMin
);
28
push_back(
"MET"
,
metMin
);
29
set(
"Muon Pt"
); set(
"MET"
);
30
wMuon_
= 0;
met_
= 0;
31
if
(
params
.exists(
"cutsToIgnore"
) ){
32
setIgnoredCuts(
params
.getParameter<std::vector<std::string> >(
"cutsToIgnore"
) );
33
}
34
retInternal_ = getBitTemplate();
35
}
37
virtual
~WSelector
() {}
39
pat::Muon
const
&
wMuon
()
const
{
return
*
wMuon_
;}
41
pat::MET
const
&
met
()
const
{
return
*
met_
; }
42
44
virtual
bool
operator()
(
edm::EventBase
const
&
event
,
pat::strbitset
&
ret
){
45
ret
.set(
false
);
46
// Handle to the muon collection
47
edm::Handle<std::vector<pat::Muon>
>
muons
;
48
// Handle to the MET collection
49
edm::Handle<std::vector<pat::MET>
>
met
;
50
// get the objects from the event
51
bool
gotMuons =
event
.getByLabel(
muonSrc_
,
muons
);
52
bool
gotMET =
event
.getByLabel(
metSrc_
,
met
);
53
// get the MET, require to be > minimum
54
if
( gotMET ){
55
met_
= &
met
->at(0);
56
if
(
met_
->
pt
() >
cut
(
metIndex_
,
double
()) || ignoreCut(
metIndex_
) )
57
passCut(
ret
,
metIndex_
);
58
}
59
// get the highest pt muon, require to have pt > minimum
60
if
( gotMuons ){
61
if
( !ignoreCut(
muonPtIndex_
) ){
62
if
(
muons
->size() > 0 ){
63
wMuon_
= &
muons
->at(0);
64
if
(
wMuon_
->
pt
() >
cut
(
muonPtIndex_
,
double
()) || ignoreCut(
muonPtIndex_
) )
65
passCut(
ret
,
muonPtIndex_
);
66
}
67
}
68
else
{
69
passCut(
ret
,
muonPtIndex_
);
70
}
71
}
72
setIgnored(
ret
);
73
return
(
bool
)
ret
;
74
}
75
76
protected
:
78
edm::InputTag
muonSrc_
;
80
edm::InputTag
metSrc_
;
82
pat::Muon
const
*
wMuon_
;
84
pat::MET
const
*
met_
;
86
index_type
muonPtIndex_
;
88
index_type
metIndex_
;
89
};
90
#endif // PhysicsTools_FWLite_WSelectorFast_h
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition:
runTheMatrix.py:355
PDWG_BPHSkim_cff.muons
muons
Definition:
PDWG_BPHSkim_cff.py:47
WSelector
Example class of an EventSelector to apply a simple W Boson selection.
Definition:
WSelector.h:16
TkAlMuonSelectors_cfi.cut
cut
Definition:
TkAlMuonSelectors_cfi.py:5
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
edm
HLT enums.
Definition:
AlignableModifier.h:19
Muon.h
WSelector::metSrc_
edm::InputTag metSrc_
met input
Definition:
WSelector.h:77
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition:
LeafCandidate.h:146
pat::Muon
Analysis-level muon class.
Definition:
Muon.h:51
edm::Handle
Definition:
AssociativeIterator.h:50
WSelector::wMuon
pat::Muon const & wMuon() const
return muon candidate of W boson
Definition:
WSelectorFast.h:39
WSelector::~WSelector
virtual ~WSelector()
destructor
Definition:
WSelectorFast.h:37
WSelector::met_
pat::MET const * met_
MET from W boson.
Definition:
WSelector.h:81
WSelector::WSelector
WSelector(edm::ParameterSet const ¶ms)
constructor
Definition:
WSelectorFast.h:21
mhtProducer_cfi.muonPtMin
muonPtMin
Definition:
mhtProducer_cfi.py:20
EventSelector.h
HLT_2018_cff.InputTag
InputTag
Definition:
HLT_2018_cff.py:79016
edm::ParameterSet
Definition:
ParameterSet.h:36
WSelector::wMuon_
pat::Muon const * wMuon_
muon candidate from W boson
Definition:
WSelector.h:79
wplusjetsAnalysis_cfi.metMin
metMin
Definition:
wplusjetsAnalysis_cfi.py:89
pat::MET
Analysis-level MET class.
Definition:
MET.h:40
WSelector::met
pat::MET const & met() const
return MET of W boson
Definition:
WSelectorFast.h:41
MET.h
WSelector::muonSrc_
edm::InputTag muonSrc_
muon input
Definition:
WSelector.h:75
pat::strbitset
Definition:
strbitset.h:23
WSelector::operator()
virtual bool operator()(edm::EventBase const &event, pat::strbitset &ret)
here is where the selection occurs
Definition:
WSelectorFast.h:44
edm::EventBase
Definition:
EventBase.h:46
EventSelector
A selector of events.
Definition:
EventSelector.h:16
event
Definition:
event.py:1
WSelector::metIndex_
index_type metIndex_
index for MET cut
Definition:
WSelectorFast.h:88
edm::InputTag
Definition:
InputTag.h:15
WSelector::muonPtIndex_
index_type muonPtIndex_
index for muon Pt cut
Definition:
WSelectorFast.h:86
Generated for CMSSW Reference Manual by
1.8.16