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 <vector>
47 #include <utility>
48 #include <tuple>
49 #include <string>
50 
51 template <typename T>
53 public:
54  explicit L1TJetsMatching(const edm::ParameterSet&);
55  ~L1TJetsMatching() override = default;
56  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58 
59 private:
60  std::pair<std::vector<T>, std::vector<T>> categorise(
61  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const;
62  std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> categoriseVBFPlus2CentralJets(
63  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const;
67  const double pt1Min_;
68  const double pt2Min_;
69  const double pt3Min_;
70  const double mjjMin_;
71  const double matchingR2_;
72 };
73 //
74 // class declaration
75 //
76 template <typename T>
77 std::pair<std::vector<T>, std::vector<T>> L1TJetsMatching<T>::categorise(
78  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const {
79  std::pair<std::vector<T>, std::vector<T>> output;
80  unsigned int i1 = 0;
81  unsigned int i2 = 0;
82  double m2jj = 0;
83  if (pfMatchedJets.size() > 1) {
84  for (unsigned int i = 0; i < pfMatchedJets.size() - 1; i++) {
85  const T& myJet1 = (pfMatchedJets)[i];
86 
87  for (unsigned int j = i + 1; j < pfMatchedJets.size(); j++) {
88  const T& myJet2 = (pfMatchedJets)[j];
89 
90  const double m2jj_test = (myJet1.p4() + myJet2.p4()).M2();
91 
92  if (m2jj_test > m2jj) {
93  m2jj = m2jj_test;
94  i1 = i;
95  i2 = j;
96  }
97  }
98  }
99 
100  const T& myJet1 = (pfMatchedJets)[i1];
101  const T& myJet2 = (pfMatchedJets)[i2];
102  const double M2jj = (Mjj >= 0. ? Mjj * Mjj : -1.);
103 
104  if ((m2jj > M2jj) && (myJet1.pt() >= pt1) && (myJet2.pt() > pt2)) {
105  output.first.push_back(myJet1);
106  output.first.push_back(myJet2);
107  }
108 
109  if ((m2jj > M2jj) && (myJet1.pt() < pt3) && (myJet1.pt() > pt2) && (myJet2.pt() > pt2)) {
110  const T& myJetTest = (pfMatchedJets)[0];
111  if (myJetTest.pt() > pt3) {
112  output.second.push_back(myJet1);
113  output.second.push_back(myJet2);
114  output.second.push_back(myJetTest);
115  }
116  }
117  }
118 
119  return output;
120 }
121 template <typename T>
122 std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> L1TJetsMatching<T>::categoriseVBFPlus2CentralJets(
123  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const { //60, 30, 50, 500
124  std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> output;
125  unsigned int i1 = 0;
126  unsigned int i2 = 0;
127 
128  double m2jj = 0;
129  if (pfMatchedJets.size() > 1) {
130  for (unsigned int i = 0; i < pfMatchedJets.size() - 1; i++) {
131  const T& myJet1 = (pfMatchedJets)[i];
132 
133  for (unsigned int j = i + 1; j < pfMatchedJets.size(); j++) {
134  const T& myJet2 = (pfMatchedJets)[j];
135 
136  const double m2jj_test = (myJet1.p4() + myJet2.p4()).M2();
137 
138  if (m2jj_test > m2jj) {
139  m2jj = m2jj_test;
140  i1 = i;
141  i2 = j;
142  }
143  }
144  }
145 
146  const T& myJet1 = (pfMatchedJets)[i1];
147  const T& myJet2 = (pfMatchedJets)[i2];
148  const double M2jj = (Mjj >= 0. ? Mjj * Mjj : -1.);
149 
150  std::vector<T> vec4jets;
151  vec4jets.reserve(4);
152  std::vector<T> vec5jets;
153  vec5jets.reserve(5);
154  std::vector<T> vec6jets;
155  vec6jets.reserve(6);
156  if (pfMatchedJets.size() > 3) {
157  if ((m2jj > M2jj) && (myJet1.pt() >= pt3) && (myJet2.pt() > pt2)) {
158  vec4jets.push_back(myJet1);
159  vec4jets.push_back(myJet2);
160 
161  for (unsigned int i = 0; i < pfMatchedJets.size(); i++) {
162  if (vec4jets.size() > 3)
163  break;
164  if (i == i1 or i == i2)
165  continue;
166  vec4jets.push_back(pfMatchedJets[i]);
167  }
168  }
169 
170  if ((m2jj > M2jj) && (myJet1.pt() < pt1) && (myJet1.pt() < pt3) && (myJet1.pt() > pt2) &&
171  (myJet2.pt() > pt2)) { //60, 30, 50, 500
172 
173  std::vector<unsigned int> idx_jets;
174  idx_jets.reserve(pfMatchedJets.size() - 2);
175 
176  for (unsigned int i = 0; i < pfMatchedJets.size(); i++) {
177  if (i == i1 || i == i2)
178  continue;
179  if (pfMatchedJets[i].pt() > pt2) {
180  idx_jets.push_back(i);
181  }
182  }
183  if (idx_jets.size() == 3) {
184  vec5jets.push_back(myJet1);
185  vec5jets.push_back(myJet2);
186  vec5jets.push_back(pfMatchedJets[idx_jets[0]]);
187  vec5jets.push_back(pfMatchedJets[idx_jets[1]]);
188  vec5jets.push_back(pfMatchedJets[idx_jets[2]]);
189 
190  } else if (idx_jets.size() > 3) {
191  vec6jets.push_back(myJet1);
192  vec6jets.push_back(myJet2);
193  vec6jets.push_back(pfMatchedJets[idx_jets[0]]);
194  vec6jets.push_back(pfMatchedJets[idx_jets[1]]);
195  vec6jets.push_back(pfMatchedJets[idx_jets[2]]);
196  vec6jets.push_back(pfMatchedJets[idx_jets[3]]);
197  }
198  }
199  }
200 
201  output = std::make_tuple(vec4jets, vec5jets, vec6jets);
202  }
203 
204  return output;
205 }
206 
207 template <typename T>
209  : jetSrc_(consumes<std::vector<T>>(iConfig.getParameter<edm::InputTag>("JetSrc"))),
210  jetTrigger_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("L1JetTrigger"))),
211  matchingMode_(iConfig.getParameter<std::string>("matchingMode")),
212  pt1Min_(iConfig.getParameter<double>("pt1Min")),
213  pt2Min_(iConfig.getParameter<double>("pt2Min")),
214  pt3Min_(iConfig.getParameter<double>("pt3Min")),
215  mjjMin_(iConfig.getParameter<double>("mjjMin")),
216  matchingR2_(iConfig.getParameter<double>("matchingR") * iConfig.getParameter<double>("matchingR")) {
217  if (matchingMode_ == "VBF") { // Default
218  produces<std::vector<T>>("TwoJets");
219  produces<std::vector<T>>("ThreeJets");
220  } else if (matchingMode_ == "VBFPlus2CentralJets") {
221  produces<std::vector<T>>("FourJets");
222  produces<std::vector<T>>("FiveJets");
223  produces<std::vector<T>>("SixJets");
224  } else {
225  throw cms::Exception("InvalidConfiguration") << "invalid value for \"matchingMode\": " << matchingMode_
226  << " (valid values are \"VBF\" and \"VBFPlus2CentralJets\")";
227  }
228 }
229 
230 template <typename T>
232  unique_ptr<std::vector<T>> pfMatchedJets(new std::vector<T>);
233 
234  // Getting HLT jets to be matched
236  iEvent.getByToken(jetSrc_, pfJets);
237 
239  iEvent.getByToken(jetTrigger_, l1TriggeredJets);
240 
241  //l1t::TauVectorRef jetCandRefVec;
242  l1t::JetVectorRef jetCandRefVec;
243  l1TriggeredJets->getObjects(trigger::TriggerL1Jet, jetCandRefVec);
244 
245  math::XYZPoint a(0., 0., 0.);
246 
247  //std::cout<<"PFsize= "<<pfJets->size()<<endl<<" L1size= "<<jetCandRefVec.size()<<std::endl;
248  for (unsigned int iJet = 0; iJet < pfJets->size(); iJet++) {
249  const T& myJet = (*pfJets)[iJet];
250  for (unsigned int iL1Jet = 0; iL1Jet < jetCandRefVec.size(); iL1Jet++) {
251  // Find the relative L2pfJets, to see if it has been reconstructed
252  // if ((iJet<3) && (iL1Jet==0)) std::cout<<myJet.p4().Pt()<<" ";
253  if ((reco::deltaR2(myJet.p4(), jetCandRefVec[iL1Jet]->p4()) < matchingR2_) && (myJet.pt() > pt2Min_)) {
254  pfMatchedJets->push_back(myJet);
255  break;
256  }
257  }
258  }
259  // order pfMatchedJets by pT
260  std::sort(pfMatchedJets->begin(), pfMatchedJets->end(), [](const T& j1, const T& j2) { return j1.pt() > j2.pt(); });
261 
262  if (matchingMode_ == "VBF") { // Default
263  std::pair<std::vector<T>, std::vector<T>> output = categorise(*pfMatchedJets, pt1Min_, pt2Min_, pt3Min_, mjjMin_);
264  auto output1 = std::make_unique<std::vector<T>>(output.first);
265  auto output2 = std::make_unique<std::vector<T>>(output.second);
266 
267  iEvent.put(std::move(output1), "TwoJets");
268  iEvent.put(std::move(output2), "ThreeJets");
269 
270  } else if (matchingMode_ == "VBFPlus2CentralJets") {
271  std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> output =
272  categoriseVBFPlus2CentralJets(*pfMatchedJets, pt1Min_, pt2Min_, pt3Min_, mjjMin_);
273  auto output1 = std::make_unique<std::vector<T>>(std::get<0>(output));
274  auto output2 = std::make_unique<std::vector<T>>(std::get<1>(output));
275  auto output3 = std::make_unique<std::vector<T>>(std::get<2>(output));
276 
277  iEvent.put(std::move(output1), "FourJets");
278  iEvent.put(std::move(output2), "FiveJets");
279  iEvent.put(std::move(output3), "SixJets");
280  }
281 }
282 
283 template <typename T>
286  desc.add<edm::InputTag>("L1JetTrigger", edm::InputTag("hltL1DiJetVBF"))->setComment("Name of trigger filter");
287  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltAK4PFJetsTightIDCorrected"))
288  ->setComment("Input collection of PFJets");
289  desc.add<std::string>("matchingMode", "VBF")
290  ->setComment("Switch from Di/tri-jet (VBF) to Multi-jet (VBFPlus2CentralJets) matching");
291  desc.add<double>("pt1Min", 110.0)->setComment("Minimal pT1 of PFJets to match");
292  desc.add<double>("pt2Min", 35.0)->setComment("Minimal pT2 of PFJets to match");
293  desc.add<double>("pt3Min", 110.0)->setComment("Minimum pT3 of PFJets to match");
294  desc.add<double>("mjjMin", 650.0)->setComment("Minimal mjj of matched PFjets");
295  desc.add<double>("matchingR", 0.5)->setComment("dR value used for matching");
296  descriptions.setComment(
297  "This module produces collection of PFJets matched to L1 Taus / Jets passing a HLT filter (Only p4 and vertex "
298  "of returned PFJets are set).");
299  descriptions.add(defaultModuleLabel<L1TJetsMatching<T>>(), desc);
300 }
301 
302 #endif
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
const double mjjMin_
L1TJetsMatching(const edm::ParameterSet &)
std::pair< std::vector< T >, std::vector< T > > categorise(const std::vector< T > &pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const
std::string defaultModuleLabel()
const edm::EDGetTokenT< std::vector< T > > jetSrc_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
~L1TJetsMatching() override=default
const double pt3Min_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::tuple< std::vector< T >, std::vector< T >, std::vector< T > > categoriseVBFPlus2CentralJets(const std::vector< T > &pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > jetTrigger_
const double matchingR2_
void setComment(std::string const &value)
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
const double pt1Min_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
double a
Definition: hdecay.h:121
std::vector< JetRef > JetVectorRef
Definition: Jet.h:14
const std::string matchingMode_
long double T
const double pt2Min_
def move(src, dest)
Definition: eostools.py:511