CMS 3D CMS Logo

PFTauL1TJetsMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
7 
8 //
9 // class declaration
10 //
12  : tauSrc_(consumes<reco::PFTauCollection>(iConfig.getParameter<edm::InputTag>("TauSrc"))),
13  L1JetSrc_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("L1JetSrc"))),
14  matchingR2_(iConfig.getParameter<double>("MatchingdR") * iConfig.getParameter<double>("MatchingdR")),
15  minTauPt_(iConfig.getParameter<double>("MinTauPt")),
16  minL1TPt_(iConfig.getParameter<double>("MinL1TPt")) {
17  produces<reco::PFTauCollection>();
18 }
20 
22  std::unique_ptr<reco::PFTauCollection> L1TmatchedPFTau(new reco::PFTauCollection);
23 
25  iEvent.getByToken(tauSrc_, taus);
26 
28  iEvent.getByToken(L1JetSrc_, L1Jets);
29 
30  l1t::JetVectorRef jetCandRefVec;
31  L1Jets->getObjects(trigger::TriggerL1Jet, jetCandRefVec);
32 
33  /* Loop over taus that must pass a certain minTauPt_ cut */
34  /* then loop over L1T jets that must pass minL1TPt_ */
35  /* and check whether they match, if yes -> include the taus in */
36  /* the new L1T matched PFTau collection */
37 
38  for (unsigned int iTau = 0; iTau < taus->size(); iTau++) {
39  bool isMatched = false;
40  if ((*taus)[iTau].pt() > minTauPt_) {
41  for (unsigned int iJet = 0; iJet < jetCandRefVec.size(); iJet++) {
42  if (jetCandRefVec[iJet]->pt() > minL1TPt_) {
43  if (reco::deltaR2((*taus)[iTau].p4(), jetCandRefVec[iJet]->p4()) < matchingR2_) {
44  isMatched = true;
45  break;
46  }
47  }
48  }
49  }
50  if (isMatched)
51  L1TmatchedPFTau->push_back((*taus)[iTau]);
52  }
53  iEvent.put(std::move(L1TmatchedPFTau));
54 }
55 
58  desc.add<edm::InputTag>("L1JetSrc", edm::InputTag("hltL1VBFDiJetOR"))
59  ->setComment("Input filter objects passing L1 seed");
60  desc.add<edm::InputTag>("TauSrc", edm::InputTag("hltSelectedPFTausTrackFindingLooseChargedIsolationAgainstMuon"))
61  ->setComment("Input collection of PFTaus");
62  desc.add<double>("MatchingdR", 0.5)->setComment("Maximum dR for matching between PFTaus and L1 filter jets");
63  desc.add<double>("MinTauPt", 20.0)->setComment("PFTaus above this pt will be considered");
64  desc.add<double>("MinL1TPt", 115.0)->setComment("L1T Objects above this pt will be considered");
65  descriptions.setComment(
66  "This module produces a collection of PFTaus matched to the leading jet passing the L1 seed filter.");
67  descriptions.add("PFTauL1TJetsMatching", desc);
68 }
edm::StreamID
Definition: StreamID.h:30
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
Tau3MuMonitor_cff.taus
taus
Definition: Tau3MuMonitor_cff.py:7
l1t::JetVectorRef
std::vector< JetRef > JetVectorRef
Definition: Jet.h:14
edm
HLT enums.
Definition: AlignableModifier.h:19
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
PFTauL1TJetsMatching::minL1TPt_
const double minL1TPt_
Definition: PFTauL1TJetsMatching.h:26
TriggerTypeDefs.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
deltaR.h
EDMException.h
PFTauL1TJetsMatching::minTauPt_
const double minTauPt_
Definition: PFTauL1TJetsMatching.h:25
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
l1t::PFTauCollection
std::vector< l1t::PFTau > PFTauCollection
Definition: PFTau.h:61
reco::PFTauCollection
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
PFTauL1TJetsMatching::~PFTauL1TJetsMatching
~PFTauL1TJetsMatching() override
Definition: PFTauL1TJetsMatching.cc:19
PFTauL1TJetsMatching.h
PFTauL1TJetsMatching::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PFTauL1TJetsMatching.cc:56
trackerHitRTTI::isMatched
bool isMatched(TrackingRecHit const &hit)
Definition: trackerHitRTTI.h:32
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
edm::ConfigurationDescriptions::setComment
void setComment(std::string const &value)
Definition: ConfigurationDescriptions.cc:48
PFTauL1TJetsMatching::PFTauL1TJetsMatching
PFTauL1TJetsMatching(const edm::ParameterSet &)
Definition: PFTauL1TJetsMatching.cc:11
HLTMonTau_cfi.L1Jets
L1Jets
Definition: HLTMonTau_cfi.py:42
iEvent
int iEvent
Definition: GenABIO.cc:224
trigger::TriggerL1Jet
Definition: TriggerTypeDefs.h:48
p4
double p4[4]
Definition: TauolaWrapper.h:92
edm::EventSetup
Definition: EventSetup.h:57
PFTauL1TJetsMatching::tauSrc_
const edm::EDGetTokenT< reco::PFTauCollection > tauSrc_
Definition: PFTauL1TJetsMatching.h:22
eostools.move
def move(src, dest)
Definition: eostools.py:511
PFTau.h
trigger
Definition: HLTPrescaleTableCond.h:8
PFTauL1TJetsMatching::L1JetSrc_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > L1JetSrc_
Definition: PFTauL1TJetsMatching.h:23
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
edm::Event
Definition: Event.h:73
PFTauL1TJetsMatching::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: PFTauL1TJetsMatching.cc:21
edm::InputTag
Definition: InputTag.h:15
PFTauL1TJetsMatching::matchingR2_
const double matchingR2_
Definition: PFTauL1TJetsMatching.h:24