CMS 3D CMS Logo

L1HLTJetsMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
10 //
11 // class decleration
12 //
13 using namespace reco;
14 using namespace std;
15 using namespace edm;
16 using namespace l1extra;
17 
19  jetSrc = consumes<edm::View<reco::Candidate> >(iConfig.getParameter<InputTag>("JetSrc"));
20  tauTrigger = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"));
21  mEt_Min = iConfig.getParameter<double>("EtMin");
22 
23  produces<CaloJetCollection>();
24 }
25 
27  using namespace edm;
28  using namespace std;
29  using namespace reco;
30  using namespace trigger;
31  using namespace l1extra;
32 
33  typedef std::vector<LeafCandidate> LeafCandidateCollection;
34 
35  unique_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
36 
37  constexpr double matchingR2 = 0.5 * 0.5;
38 
39  //Getting HLT jets to be matched
41  iEvent.getByToken(jetSrc, tauJets);
42 
43  // std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
44 
46  iEvent.getByToken(tauTrigger, l1TriggeredTaus);
47 
48  std::vector<l1extra::L1JetParticleRef> tauCandRefVec;
49  std::vector<l1extra::L1JetParticleRef> jetCandRefVec;
50 
51  l1TriggeredTaus->getObjects(trigger::TriggerL1TauJet, tauCandRefVec);
52  l1TriggeredTaus->getObjects(trigger::TriggerL1CenJet, jetCandRefVec);
53  math::XYZPoint a(0., 0., 0.);
55 
56  for (unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++) {
57  for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
58  //Find the relative L2TauJets, to see if it has been reconstructed
59  const Candidate& myJet = (*tauJets)[iJet];
60  double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
61  if (deltaR2 < matchingR2) {
62  // LeafCandidate myLC(myJet);
63  CaloJet myCaloJet(myJet.p4(), a, f);
64  if (myJet.pt() > mEt_Min) {
65  // tauL2LC->push_back(myLC);
66  tauL2jets->push_back(myCaloJet);
67  }
68  break;
69  }
70  }
71  }
72 
73  for (unsigned int iL1Tau = 0; iL1Tau < jetCandRefVec.size(); iL1Tau++) {
74  for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
75  const Candidate& myJet = (*tauJets)[iJet];
76  //Find the relative L2TauJets, to see if it has been reconstructed
77  double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
78  if (deltaR2 < matchingR2) {
79  // LeafCandidate myLC(myJet);
80  CaloJet myCaloJet(myJet.p4(), a, f);
81  if (myJet.pt() > mEt_Min) {
82  //tauL2LC->push_back(myLC);
83  tauL2jets->push_back(myCaloJet);
84  }
85  break;
86  }
87  }
88  }
89 
90  //std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;
91 
92  iEvent.put(std::move(tauL2jets));
93  // iEvent.put(std::move(tauL2LC));
94 }
95 
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Jets made from CaloTowers.
Definition: CaloJet.h:27
virtual double pt() const =0
transverse momentum
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int iEvent
Definition: GenABIO.cc:224
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
L1HLTJetsMatching(const edm::ParameterSet &)
def move(src, dest)
Definition: eostools.py:511
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector