CMS 3D CMS Logo

L1HLTTauMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
8 
9 //
10 // class decleration
11 //
12 using namespace reco;
13 using namespace std;
14 using namespace edm;
15 using namespace l1extra;
16 
18  : jetSrc(consumes<PFTauCollection>(iConfig.getParameter<InputTag>("JetSrc"))),
19  tauTrigger(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"))),
20  mEt_Min(iConfig.getParameter<double>("EtMin")) {
21  produces<PFTauCollection>();
22 }
24 
26  using namespace edm;
27  using namespace std;
28  using namespace reco;
29  using namespace trigger;
30  using namespace l1extra;
31 
32  unique_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
33 
34  double deltaR = 1.0;
35  double matchingR = 0.5;
36  //Getting HLT jets to be matched
38  iEvent.getByToken(jetSrc, tauJets);
39 
40  // std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
41 
43  iEvent.getByToken(tauTrigger, l1TriggeredTaus);
44 
45  vector<l1extra::L1JetParticleRef> tauCandRefVec;
46  vector<l1extra::L1JetParticleRef> jetCandRefVec;
47  l1TriggeredTaus->getObjects(trigger::TriggerL1TauJet, tauCandRefVec);
48  l1TriggeredTaus->getObjects(trigger::TriggerL1CenJet, jetCandRefVec);
49 
50  math::XYZPoint a(0., 0., 0.);
51 
52  for (unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++) {
53  for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
54  //Find the relative L2TauJets, to see if it has been reconstructed
55  const PFTau& myJet = (*tauJets)[iJet];
56  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
57  if (deltaR < matchingR) {
58  // LeafCandidate myLC(myJet);
59  if (myJet.leadChargedHadrCand().isNonnull()) {
60  a = myJet.leadChargedHadrCand()->vertex();
61  }
62  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
63  if (myJet.pt() > mEt_Min) {
64  // tauL2LC->push_back(myLC);
65  tauL2jets->push_back(myPFTau);
66  }
67  break;
68  }
69  }
70  }
71 
72  for (unsigned int iL1Tau = 0; iL1Tau < jetCandRefVec.size(); iL1Tau++) {
73  for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
74  const PFTau& myJet = (*tauJets)[iJet];
75  //Find the relative L2TauJets, to see if it has been reconstructed
76  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
77  if (deltaR < matchingR) {
78  // LeafCandidate myLC(myJet);
79  if (myJet.leadChargedHadrCand().isNonnull()) {
80  a = myJet.leadChargedHadrCand()->vertex();
81  }
82 
83  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
84  if (myJet.pt() > mEt_Min) {
85  // tauL2LC->push_back(myLC);
86  tauL2jets->push_back(myPFTau);
87  }
88  break;
89  }
90  }
91  }
92 
93  //std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;
94 
95  iEvent.put(std::move(tauL2jets));
96  // iEvent.put(std::move(tauL2LC));
97 }
98 
101  desc.add<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleIsoTau40er"))
102  ->setComment("Name of trigger filter");
103  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltSelectedPFTausTrackPt1MediumIsolationReg"))
104  ->setComment("Input collection of PFTaus");
105  desc.add<double>("EtMin", 0.0)->setComment("Minimal pT of PFTau to match");
106  descriptions.setComment(
107  "This module produces collection of PFTaus matched to L1ExtraTaus/Jets passing a HLT filter (Only p4 and vertex "
108  "of returned PFTaus are set).");
109  descriptions.add("L1HLTJetsMatching", desc);
110 }
edm::StreamID
Definition: StreamID.h:30
L1HLTTauMatching::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1HLTTauMatching.cc:99
L1HLTTauMatching::tauTrigger
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauTrigger
Definition: L1HLTTauMatching.h:30
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
edm
HLT enums.
Definition: AlignableModifier.h:19
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:45
HLT_2018_cff.matchingR
matchingR
Definition: HLT_2018_cff.py:66730
edm::Handle
Definition: AssociativeIterator.h:50
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:593
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
L1HLTTauMatching::mEt_Min
const double mEt_Min
Definition: L1HLTTauMatching.h:31
L1HLTTauMatching::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: L1HLTTauMatching.cc:25
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
reco::PFTau::leadChargedHadrCand
const CandidatePtr & leadChargedHadrCand() const
Definition: PFTau.cc:62
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
L1JetParticleFwd.h
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
a
double a
Definition: hdecay.h:119
edm::ConfigurationDescriptions::setComment
void setComment(std::string const &value)
Definition: ConfigurationDescriptions.cc:48
reco::Candidate::vertex
virtual const Point & vertex() const =0
vertex position
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::LeafCandidate::p4
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:114
L1HLTTauMatching::~L1HLTTauMatching
~L1HLTTauMatching() override
Definition: L1HLTTauMatching.cc:23
edm::EventSetup
Definition: EventSetup.h:57
electronAnalyzer_cfi.DeltaR
DeltaR
Definition: electronAnalyzer_cfi.py:33
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
trigger::TriggerL1TauJet
Definition: TriggerTypeDefs.h:35
L1HLTTauMatching::jetSrc
const edm::EDGetTokenT< reco::PFTauCollection > jetSrc
Definition: L1HLTTauMatching.h:29
PFTau.h
edm::Ptr::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
L1HLTTauMatching::L1HLTTauMatching
L1HLTTauMatching(const edm::ParameterSet &)
Definition: L1HLTTauMatching.cc:17
trigger
Definition: HLTPrescaleTableCond.h:8
trigger::TriggerL1CenJet
Definition: TriggerTypeDefs.h:33
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
l1extra
Definition: L1EmParticle.h:26
edm::Event
Definition: Event.h:73
L1JetParticle.h
DisplacedJet_Monitor_cff.jetSrc
jetSrc
Definition: DisplacedJet_Monitor_cff.py:79
edm::InputTag
Definition: InputTag.h:15
L1HLTTauMatching.h