CMS 3D CMS Logo

L1THLTTauMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
6 
7 //
8 // class decleration
9 //
10 using namespace reco;
11 using namespace std;
12 using namespace edm;
13 using namespace trigger;
14 
16  : jetSrc(consumes<PFTauCollection>(iConfig.getParameter<InputTag>("JetSrc"))),
17  tauTrigger(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"))),
18  mEt_Min(iConfig.getParameter<double>("EtMin")) {
19  produces<PFTauCollection>();
20 }
22 
24  unique_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
25 
26  double deltaR = 1.0;
27  double matchingR = 0.5;
28 
29  // Getting HLT jets to be matched
31  iEvent.getByToken(jetSrc, tauJets);
32 
34  iEvent.getByToken(tauTrigger, l1TriggeredTaus);
35 
36  l1t::TauVectorRef tauCandRefVec;
37  l1TriggeredTaus->getObjects(trigger::TriggerL1Tau, tauCandRefVec);
38 
39  math::XYZPoint a(0., 0., 0.);
40 
41  for (unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++) {
42  for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
43  // Find the relative L2TauJets, to see if it has been reconstructed
44  const PFTau& myJet = (*tauJets)[iJet];
45  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
46  if (deltaR < matchingR) {
47  if (myJet.leadChargedHadrCand().isNonnull()) {
48  a = myJet.leadChargedHadrCand()->vertex();
49  }
50  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
51  if (myJet.pt() > mEt_Min) {
52  tauL2jets->push_back(myPFTau);
53  }
54  break;
55  }
56  }
57  }
58 
59  iEvent.put(std::move(tauL2jets));
60 }
61 
64  desc.add<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleIsoTau40er"))
65  ->setComment("Name of trigger filter");
66  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltSelectedPFTausTrackPt1MediumIsolationReg"))
67  ->setComment("Input collection of PFTaus");
68  desc.add<double>("EtMin", 0.0)->setComment("Minimal pT of PFTau to match");
69  descriptions.setComment(
70  "This module produces collection of PFTaus matched to L1 Taus / Jets passing a HLT filter (Only p4 and vertex of "
71  "returned PFTaus are set).");
72  descriptions.add("L1THLTTauMatching", desc);
73 }
edm::StreamID
Definition: StreamID.h:30
trigger::TriggerFilterObjectWithRefs
Definition: TriggerFilterObjectWithRefs.h:35
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
TriggerTypeDefs.h
reco::PFTau
Definition: PFTau.h:36
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition: LeafCandidate.h:146
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
edm::Handle
Definition: AssociativeIterator.h:50
L1THLTTauMatching::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1THLTTauMatching.cc:62
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:590
EDMException.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
reco::PFTauCollection
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
reco::PFTau::leadChargedHadrCand
const CandidatePtr & leadChargedHadrCand() const
Definition: PFTau.cc:63
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
edm::ParameterSet
Definition: ParameterSet.h:47
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
a
double a
Definition: hdecay.h:119
L1THLTTauMatching.h
trigger::TriggerL1Tau
Definition: TriggerTypeDefs.h:49
L1THLTTauMatching::mEt_Min
const double mEt_Min
Definition: L1THLTTauMatching.h:30
edm::ConfigurationDescriptions::setComment
void setComment(std::string const &value)
Definition: ConfigurationDescriptions.cc:48
L1THLTTauMatching::tauTrigger
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauTrigger
Definition: L1THLTTauMatching.h:29
reco::Candidate::vertex
virtual const Point & vertex() const =0
vertex position
L1THLTTauMatching::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: L1THLTTauMatching.cc:23
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::LeafCandidate::p4
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:114
edm::EventSetup
Definition: EventSetup.h:58
electronAnalyzer_cfi.DeltaR
DeltaR
Definition: electronAnalyzer_cfi.py:33
HLT_FULL_cff.matchingR
matchingR
Definition: HLT_FULL_cff.py:68718
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
PFTau.h
l1t::TauVectorRef
std::vector< TauRef > TauVectorRef
Definition: Tau.h:14
edm::Ptr::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
L1THLTTauMatching::~L1THLTTauMatching
~L1THLTTauMatching() override
Definition: L1THLTTauMatching.cc:21
trigger
Definition: HLTPrescaleTableCond.h:8
edm::Event
Definition: Event.h:73
DisplacedJet_Monitor_cff.jetSrc
jetSrc
Definition: DisplacedJet_Monitor_cff.py:79
edm::InputTag
Definition: InputTag.h:15
L1THLTTauMatching::jetSrc
const edm::EDGetTokenT< reco::PFTauCollection > jetSrc
Definition: L1THLTTauMatching.h:28
L1THLTTauMatching::L1THLTTauMatching
L1THLTTauMatching(const edm::ParameterSet &)
Definition: L1THLTTauMatching.cc:15