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;
56  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58  std::pair<std::vector<T>, std::vector<T>> categorise(
59  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const;
60  std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> categoriseVBFPlus2CentralJets(
61  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const;
62 
63 private:
67  const double pt1Min_;
68  const double pt2Min_;
69  const double pt3Min_;
70  const double mjjMin_;
71  const double matchingR_;
72  const double matchingR2_;
73 };
74 //
75 // class declaration
76 //
77 template <typename T>
78 std::pair<std::vector<T>, std::vector<T>> L1TJetsMatching<T>::categorise(
79  const std::vector<T>& pfMatchedJets, double pt1, double pt2, double pt3, double Mjj) const {
80  std::pair<std::vector<T>, std::vector<T>> output;
81  unsigned int i1 = 0;
82  unsigned int i2 = 0;
83  double mjj = 0;
84  if (pfMatchedJets.size() > 1) {
85  for (unsigned int i = 0; i < pfMatchedJets.size() - 1; i++) {
86  const T& myJet1 = (pfMatchedJets)[i];
87 
88  for (unsigned int j = i + 1; j < pfMatchedJets.size(); j++) {
89  const T& myJet2 = (pfMatchedJets)[j];
90 
91  const double mjj_test = (myJet1.p4() + myJet2.p4()).M();
92 
93  if (mjj_test > mjj) {
94  mjj = mjj_test;
95  i1 = i;
96  i2 = j;
97  }
98  }
99  }
100 
101  const T& myJet1 = (pfMatchedJets)[i1];
102  const T& myJet2 = (pfMatchedJets)[i2];
103 
104  if ((mjj > Mjj) && (myJet1.pt() >= pt1) && (myJet2.pt() > pt2)) {
105  output.first.push_back(myJet1);
106  output.first.push_back(myJet2);
107  }
108 
109  if ((mjj > Mjj) && (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 mjj = 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 mjj_test = (myJet1.p4() + myJet2.p4()).M();
137 
138  if (mjj_test > mjj) {
139  mjj = mjj_test;
140  i1 = i;
141  i2 = j;
142  }
143  }
144  }
145 
146  const T& myJet1 = (pfMatchedJets)[i1];
147  const T& myJet2 = (pfMatchedJets)[i2];
148 
149  std::vector<T> vec4jets;
150  vec4jets.reserve(4);
151  std::vector<T> vec5jets;
152  vec5jets.reserve(5);
153  std::vector<T> vec6jets;
154  vec6jets.reserve(6);
155  if (pfMatchedJets.size() > 3) {
156  if ((mjj > Mjj) && (myJet1.pt() >= pt3) && (myJet2.pt() > pt2)) {
157  vec4jets.push_back(myJet1);
158  vec4jets.push_back(myJet2);
159 
160  for (unsigned int i = 0; i < pfMatchedJets.size(); i++) {
161  if (vec4jets.size() > 3)
162  break;
163  if (i == i1 or i == i2)
164  continue;
165  vec4jets.push_back(pfMatchedJets[i]);
166  }
167  }
168 
169  if ((mjj > Mjj) && (myJet1.pt() < pt1) && (myJet1.pt() < pt3) && (myJet1.pt() > pt2) &&
170  (myJet2.pt() > pt2)) { //60, 30, 50, 500
171 
172  std::vector<unsigned int> idx_jets;
173  idx_jets.reserve(pfMatchedJets.size() - 2);
174 
175  for (unsigned int i = 0; i < pfMatchedJets.size(); i++) {
176  if (i == i1 || i == i2)
177  continue;
178  if (pfMatchedJets[i].pt() > pt2) {
179  idx_jets.push_back(i);
180  }
181  }
182  if (idx_jets.size() == 3) {
183  vec5jets.push_back(myJet1);
184  vec5jets.push_back(myJet2);
185  vec5jets.push_back(pfMatchedJets[idx_jets[0]]);
186  vec5jets.push_back(pfMatchedJets[idx_jets[1]]);
187  vec5jets.push_back(pfMatchedJets[idx_jets[2]]);
188 
189  } else if (idx_jets.size() > 3) {
190  vec6jets.push_back(myJet1);
191  vec6jets.push_back(myJet2);
192  vec6jets.push_back(pfMatchedJets[idx_jets[0]]);
193  vec6jets.push_back(pfMatchedJets[idx_jets[1]]);
194  vec6jets.push_back(pfMatchedJets[idx_jets[2]]);
195  vec6jets.push_back(pfMatchedJets[idx_jets[3]]);
196  }
197  }
198  }
199 
200  output = std::make_tuple(vec4jets, vec5jets, vec6jets);
201  }
202 
203  return output;
204 }
205 
206 template <typename T>
208  : jetSrc_(consumes<std::vector<T>>(iConfig.getParameter<edm::InputTag>("JetSrc"))),
209  jetTrigger_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("L1JetTrigger"))),
210  matchingMode_(iConfig.getParameter<std::string>("matchingMode")),
211  pt1Min_(iConfig.getParameter<double>("pt1Min")),
212  pt2Min_(iConfig.getParameter<double>("pt2Min")),
213  pt3Min_(iConfig.getParameter<double>("pt3Min")),
214  mjjMin_(iConfig.getParameter<double>("mjjMin")),
215  matchingR_(iConfig.getParameter<double>("matchingR")),
216  matchingR2_(matchingR_ * 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 template <typename T>
231 
232 template <typename T>
234  unique_ptr<std::vector<T>> pfMatchedJets(new std::vector<T>);
235 
236  // Getting HLT jets to be matched
238  iEvent.getByToken(jetSrc_, pfJets);
239 
241  iEvent.getByToken(jetTrigger_, l1TriggeredJets);
242 
243  //l1t::TauVectorRef jetCandRefVec;
244  l1t::JetVectorRef jetCandRefVec;
245  l1TriggeredJets->getObjects(trigger::TriggerL1Jet, jetCandRefVec);
246 
247  math::XYZPoint a(0., 0., 0.);
248 
249  //std::cout<<"PFsize= "<<pfJets->size()<<endl<<" L1size= "<<jetCandRefVec.size()<<std::endl;
250  for (unsigned int iJet = 0; iJet < pfJets->size(); iJet++) {
251  const T& myJet = (*pfJets)[iJet];
252  for (unsigned int iL1Jet = 0; iL1Jet < jetCandRefVec.size(); iL1Jet++) {
253  // Find the relative L2pfJets, to see if it has been reconstructed
254  // if ((iJet<3) && (iL1Jet==0)) std::cout<<myJet.p4().Pt()<<" ";
255  if ((reco::deltaR2(myJet.p4(), jetCandRefVec[iL1Jet]->p4()) < matchingR2_) && (myJet.pt() > pt2Min_)) {
256  pfMatchedJets->push_back(myJet);
257  break;
258  }
259  }
260  }
261  // order pfMatchedJets by pT
262  std::sort(pfMatchedJets->begin(), pfMatchedJets->end(), [](const T& j1, const T& j2) { return j1.pt() > j2.pt(); });
263 
264  if (matchingMode_ == "VBF") { // Default
265  std::pair<std::vector<T>, std::vector<T>> output = categorise(*pfMatchedJets, pt1Min_, pt2Min_, pt3Min_, mjjMin_);
266  auto output1 = std::make_unique<std::vector<T>>(output.first);
267  auto output2 = std::make_unique<std::vector<T>>(output.second);
268 
269  iEvent.put(std::move(output1), "TwoJets");
270  iEvent.put(std::move(output2), "ThreeJets");
271 
272  } else if (matchingMode_ == "VBFPlus2CentralJets") {
273  std::tuple<std::vector<T>, std::vector<T>, std::vector<T>> output =
274  categoriseVBFPlus2CentralJets(*pfMatchedJets, pt1Min_, pt2Min_, pt3Min_, mjjMin_);
275  auto output1 = std::make_unique<std::vector<T>>(std::get<0>(output));
276  auto output2 = std::make_unique<std::vector<T>>(std::get<1>(output));
277  auto output3 = std::make_unique<std::vector<T>>(std::get<2>(output));
278 
279  iEvent.put(std::move(output1), "FourJets");
280  iEvent.put(std::move(output2), "FiveJets");
281  iEvent.put(std::move(output3), "SixJets");
282  }
283 }
284 
285 template <typename T>
288  desc.add<edm::InputTag>("L1JetTrigger", edm::InputTag("hltL1DiJetVBF"))->setComment("Name of trigger filter");
289  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltAK4PFJetsTightIDCorrected"))
290  ->setComment("Input collection of PFJets");
291  desc.add<std::string>("matchingMode", "VBF")
292  ->setComment("Switch from Di/tri-jet (VBF) to Multi-jet (VBFPlus2CentralJets) matching");
293  desc.add<double>("pt1Min", 110.0)->setComment("Minimal pT1 of PFJets to match");
294  desc.add<double>("pt2Min", 35.0)->setComment("Minimal pT2 of PFJets to match");
295  desc.add<double>("pt3Min", 110.0)->setComment("Minimum pT3 of PFJets to match");
296  desc.add<double>("mjjMin", 650.0)->setComment("Minimal mjj of matched PFjets");
297  desc.add<double>("matchingR", 0.5)->setComment("dR value used for matching");
298  descriptions.setComment(
299  "This module produces collection of PFJets matched to L1 Taus / Jets passing a HLT filter (Only p4 and vertex "
300  "of returned PFJets are set).");
301  descriptions.add(defaultModuleLabel<L1TJetsMatching<T>>(), desc);
302 }
303 
304 #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
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:119
std::vector< JetRef > JetVectorRef
Definition: Jet.h:14
const std::string matchingMode_
const double matchingR_
long double T
const double pt2Min_
def move(src, dest)
Definition: eostools.py:511