CMS 3D CMS Logo

L1TJetsMatching.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoTauTag/HLTProducers
4 // Class: L1TJetsMatching
5 //
16 //
17 // Original Author: Vukasin Milosevic
18 // Created: Thu, 01 Jun 2017 17:23:00 GMT
19 //
20 //
21 
22 #ifndef RecoTauTag_HLTProducers_L1TJetsMatching_h
23 #define RecoTauTag_HLTProducers_L1TJetsMatching_h
24 
25 // user include files
37 
38 #include "Math/GenVector/VectorUtil.h"
42 
45 
46 #include <map>
47 #include <vector>
48 
49 template <typename T>
51 public:
52  explicit L1TJetsMatching(const edm::ParameterSet&);
53  ~L1TJetsMatching() override;
54  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56 
57 private:
60  const double pt1Min_;
61  const double pt2Min_;
62  const double pt3Min_;
63  const double mjjMin_;
64  const double matchingR_;
65  const double matchingR2_;
66 };
67 //
68 // class decleration
69 //
70 template <typename T>
71 std::pair<std::vector<T>, std::vector<T>> categorise(
72  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) {
73  std::pair<std::vector<T>, std::vector<T>> output;
74  unsigned int i1 = 0;
75  unsigned int i2 = 0;
76  double mjj = 0;
77  if (pfMatchedJets.size() > 1) {
78  for (unsigned int i = 0; i < pfMatchedJets.size() - 1; i++) {
79  const T& myJet1 = (pfMatchedJets)[i];
80 
81  for (unsigned int j = i + 1; j < pfMatchedJets.size(); j++) {
82  const T& myJet2 = (pfMatchedJets)[j];
83 
84  const double mjj_test = (myJet1.p4() + myJet2.p4()).M();
85 
86  if (mjj_test > mjj) {
87  mjj = mjj_test;
88  i1 = i;
89  i2 = j;
90  }
91  }
92  }
93 
94  const T& myJet1 = (pfMatchedJets)[i1];
95  const T& myJet2 = (pfMatchedJets)[i2];
96 
97  if ((mjj > Mjj) && (myJet1.pt() >= pt1) && (myJet2.pt() > pt2)) {
98  output.first.push_back(myJet1);
99  output.first.push_back(myJet2);
100  }
101 
102  if ((mjj > Mjj) && (myJet1.pt() < pt3) && (myJet1.pt() > pt2) && (myJet2.pt() > pt2)) {
103  const T& myJetTest = (pfMatchedJets)[0];
104  if (myJetTest.pt() > pt3) {
105  output.second.push_back(myJet1);
106  output.second.push_back(myJet2);
107  output.second.push_back(myJetTest);
108  }
109  }
110  }
111 
112  return output;
113 }
114 template <typename T>
116  : jetSrc_(consumes<std::vector<T>>(iConfig.getParameter<edm::InputTag>("JetSrc"))),
117  jetTrigger_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("L1JetTrigger"))),
118  pt1Min_(iConfig.getParameter<double>("pt1Min")),
119  pt2Min_(iConfig.getParameter<double>("pt2Min")),
120  pt3Min_(iConfig.getParameter<double>("pt3Min")),
121  mjjMin_(iConfig.getParameter<double>("mjjMin")),
122  matchingR_(iConfig.getParameter<double>("matchingR")),
123  matchingR2_(matchingR_ * matchingR_) {
124  produces<std::vector<T>>("TwoJets");
125  produces<std::vector<T>>("ThreeJets");
126 }
127 template <typename T>
129 
130 template <typename T>
132  unique_ptr<std::vector<T>> pfMatchedJets(new std::vector<T>);
133  std::pair<std::vector<T>, std::vector<T>> output;
134 
135  // Getting HLT jets to be matched
137  iEvent.getByToken(jetSrc_, pfJets);
138 
140  iEvent.getByToken(jetTrigger_, l1TriggeredJets);
141 
142  //l1t::TauVectorRef jetCandRefVec;
143  l1t::JetVectorRef jetCandRefVec;
144  l1TriggeredJets->getObjects(trigger::TriggerL1Jet, jetCandRefVec);
145 
146  math::XYZPoint a(0., 0., 0.);
147 
148  //std::cout<<"PFsize= "<<pfJets->size()<<endl<<" L1size= "<<jetCandRefVec.size()<<std::endl;
149  for (unsigned int iJet = 0; iJet < pfJets->size(); iJet++) {
150  const T& myJet = (*pfJets)[iJet];
151  for (unsigned int iL1Jet = 0; iL1Jet < jetCandRefVec.size(); iL1Jet++) {
152  // Find the relative L2pfJets, to see if it has been reconstructed
153  // if ((iJet<3) && (iL1Jet==0)) std::cout<<myJet.p4().Pt()<<" ";
154  if ((reco::deltaR2(myJet.p4(), jetCandRefVec[iL1Jet]->p4()) < matchingR2_) && (myJet.pt() > pt2Min_)) {
155  pfMatchedJets->push_back(myJet);
156  break;
157  }
158  }
159  }
160 
161  output = categorise(*pfMatchedJets, pt1Min_, pt2Min_, pt3Min_, mjjMin_);
162  unique_ptr<std::vector<T>> output1(new std::vector<T>(output.first));
163  unique_ptr<std::vector<T>> output2(new std::vector<T>(output.second));
164 
165  iEvent.put(std::move(output1), "TwoJets");
166  iEvent.put(std::move(output2), "ThreeJets");
167 }
168 template <typename T>
171  desc.add<edm::InputTag>("L1JetTrigger", edm::InputTag("hltL1DiJetVBF"))->setComment("Name of trigger filter");
172  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltAK4PFJetsTightIDCorrected"))
173  ->setComment("Input collection of PFJets");
174  desc.add<double>("pt1Min", 110.0)->setComment("Minimal pT1 of PFJets to match");
175  desc.add<double>("pt2Min", 35.0)->setComment("Minimal pT2 of PFJets to match");
176  desc.add<double>("pt3Min", 110.0)->setComment("Minimum pT3 of PFJets to match");
177  desc.add<double>("mjjMin", 650.0)->setComment("Minimal mjj of matched PFjets");
178  desc.add<double>("matchingR", 0.5)->setComment("dR value used for matching");
179  descriptions.setComment(
180  "This module produces collection of PFJetss matched to L1 Taus / Jets passing a HLT filter (Only p4 and vertex "
181  "of returned PFJetss are set).");
182  descriptions.add(defaultModuleLabel<L1TJetsMatching<T>>(), desc);
183 }
184 
185 #endif
defaultModuleLabel.h
ConfigurationDescriptions.h
edm::StreamID
Definition: StreamID.h:30
testProducerWithPsetDescEmpty_cfi.i2
i2
Definition: testProducerWithPsetDescEmpty_cfi.py:46
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
PFTauFwd.h
skim900GeV_StreamA_MinBiasPD_cfg.output1
output1
Definition: skim900GeV_StreamA_MinBiasPD_cfg.py:212
L1TJetsMatching::pt3Min_
const double pt3Min_
Definition: L1TJetsMatching.h:62
L1TJetsMatching::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1TJetsMatching.h:169
L1TJetsMatching::L1TJetsMatching
L1TJetsMatching(const edm::ParameterSet &)
Definition: L1TJetsMatching.h:115
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
edm::EDGetTokenT
Definition: EDGetToken.h:33
l1t::JetVectorRef
std::vector< JetRef > JetVectorRef
Definition: Jet.h:14
edm
HLT enums.
Definition: AlignableModifier.h:19
L1TJetsMatching::jetSrc_
const edm::EDGetTokenT< std::vector< T > > jetSrc_
Definition: L1TJetsMatching.h:58
PFJet.h
testProducerWithPsetDescEmpty_cfi.i1
i1
Definition: testProducerWithPsetDescEmpty_cfi.py:45
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
PFJetCollection.h
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
TriggerTypeDefs.h
TriggerFilterObjectWithRefs.h
edm::Handle
Definition: AssociativeIterator.h:50
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:452
L1TJetsMatching::matchingR_
const double matchingR_
Definition: L1TJetsMatching.h:64
L1TJetsMatching::jetTrigger_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > jetTrigger_
Definition: L1TJetsMatching.h:59
EDMException.h
categorise
std::pair< std::vector< T >, std::vector< T > > categorise(const std::vector< T > &pfMatchedJets, double pt1, double pt2, double pt3, double Mjj)
Definition: L1TJetsMatching.h:71
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
L1TJetsMatching::pt1Min_
const double pt1Min_
Definition: L1TJetsMatching.h:60
HLT_FULL_cff.pt1
pt1
Definition: HLT_FULL_cff.py:9870
ParameterSetDescription.h
L1TJetsMatching::~L1TJetsMatching
~L1TJetsMatching() override
Definition: L1TJetsMatching.h:128
edm::global::EDProducer
Definition: EDProducer.h:32
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
defaultModuleLabel
std::string defaultModuleLabel()
Definition: defaultModuleLabel.h:16
Event.h
deltaR.h
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
L1TJetsMatching::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: L1TJetsMatching.h:131
iEvent
int iEvent
Definition: GenABIO.cc:224
trigger::TriggerL1Jet
Definition: TriggerTypeDefs.h:48
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
L1TJetsMatching::matchingR2_
const double matchingR2_
Definition: L1TJetsMatching.h:65
edm::EventSetup
Definition: EventSetup.h:57
skim900GeV_StreamA_MinBiasPD_cfg.output2
output2
Definition: skim900GeV_StreamA_MinBiasPD_cfg.py:226
HLT_FULL_cff.pt2
pt2
Definition: HLT_FULL_cff.py:9872
InputTag.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
L1TJetsMatching
Definition: L1TJetsMatching.h:50
Frameworkfwd.h
T
long double T
Definition: Basic3DVectorLD.h:48
trigger
Definition: HLTPrescaleTableCond.h:8
ParameterSet.h
EDProducer.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
L1TJetsMatching::mjjMin_
const double mjjMin_
Definition: L1TJetsMatching.h:63
L1TJetsMatching::pt2Min_
const double pt2Min_
Definition: L1TJetsMatching.h:61
edm::InputTag
Definition: InputTag.h:15
pfJetBenchmark_cfi.pfJets
pfJets
Definition: pfJetBenchmark_cfi.py:4