CMS 3D CMS Logo

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