CMS 3D CMS Logo

L1TrackSelectionProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/L1TTrackMatch
4 // Class: L1TrackSelectionProducer
5 //
18 //
19 // Original Author: Alexx Perloff
20 // Created: Thu, 16 Dec 2021 19:02:50 GMT
21 //
22 //
23 
24 // system include files
25 #include <algorithm>
26 #include <memory>
27 #include <string>
28 #include <vector>
29 
30 // Xilinx HLS includes
31 #include <ap_fixed.h>
32 #include <ap_int.h>
33 
34 // user include files
59 
60 //
61 // class declaration
62 //
63 
65 public:
67  ~L1TrackSelectionProducer() override;
68 
69  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
70 
71 private:
72  // ----------constants, enums and typedefs ---------
73  // Relevant constants for the converted track word
75  kPtSize = TTTrack_TrackWord::TrackBitWidths::kRinvSize - 1, // Width of pt
76  kPtMagSize = 9, // Width of pt magnitude (unsigned)
77  kEtaSize = TTTrack_TrackWord::TrackBitWidths::kTanlSize, // Width of eta
78  kEtaMagSize = 3, // Width of eta magnitude (signed)
79  };
80 
82  typedef std::vector<L1Track> TTTrackCollection;
85  typedef std::vector<TTTrackRef> TTTrackRefCollection;
86  typedef std::unique_ptr<TTTrackRefCollection> TTTrackRefCollectionUPtr;
87 
88  // ----------member functions ----------------------
89  void printDebugInfo(const TTTrackCollectionHandle& l1TracksHandle,
90  const TTTrackRefCollectionUPtr& vTTTrackOutput,
91  const TTTrackRefCollectionUPtr& vTTTrackEmulationOutput,
92  const TTTrackRefCollectionUPtr& vTTTrackAssociatedOutput,
93  const TTTrackRefCollectionUPtr& vTTTrackAssociatedEmulationOutput) const;
94  void printTrackInfo(edm::LogInfo& log, const L1Track& track, bool printEmulation = false) const;
95  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
96 
97  // ----------selectors -----------------------------
98  // Based on recommendations from https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideGenericSelectors
101  TTTrackPtMinSelector(const edm::ParameterSet& cfg) : ptMin_(cfg.template getParameter<double>("ptMin")) {}
102  bool operator()(const L1Track& t) const { return t.momentum().perp() >= ptMin_; }
103 
104  private:
105  double ptMin_;
106  };
109  TTTrackWordPtMinSelector(const edm::ParameterSet& cfg) : ptMin_(cfg.template getParameter<double>("ptMin")) {}
110  bool operator()(const L1Track& t) const {
111  ap_uint<TrackBitWidths::kPtSize> ptEmulationBits = t.getTrackWord()(
112  TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1, TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
113  ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize> ptEmulation;
114  ptEmulation.V = ptEmulationBits.range();
115  //edm::LogInfo("L1TrackSelectionProducer") << "produce::Emulation track properties::ap_uint(bits) = " << ptEmulationBits.to_string(2)
116  // << " ap_ufixed(bits) = " << ptEmulation.to_string(2) << " ap_ufixed(float) = " << ptEmulation.to_double();
117  return ptEmulation.to_double() >= ptMin_;
118  }
119 
120  private:
121  double ptMin_;
122  };
126  : absEtaMax_(cfg.template getParameter<double>("absEtaMax")) {}
127  bool operator()(const L1Track& t) const { return std::abs(t.momentum().eta()) <= absEtaMax_; }
128 
129  private:
130  double absEtaMax_;
131  };
135  : absEtaMax_(cfg.template getParameter<double>("absEtaMax")) {}
136  bool operator()(const L1Track& t) const {
137  TTTrack_TrackWord::tanl_t etaEmulationBits = t.getTanlWord();
138  ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize> etaEmulation;
139  etaEmulation.V = etaEmulationBits.range();
140  return std::abs(etaEmulation.to_double()) <= absEtaMax_;
141  }
142 
143  private:
144  double absEtaMax_;
145  };
148  TTTrackAbsZ0MaxSelector(const edm::ParameterSet& cfg) : absZ0Max_(cfg.template getParameter<double>("absZ0Max")) {}
149  bool operator()(const L1Track& t) const { return std::abs(t.z0()) <= absZ0Max_; }
150 
151  private:
152  double absZ0Max_;
153  };
157  : absZ0Max_(cfg.template getParameter<double>("absZ0Max")) {}
158  bool operator()(const L1Track& t) const { return std::abs(t.getZ0()) <= absZ0Max_; }
159 
160  private:
161  double absZ0Max_;
162  };
166  : nStubsMin_(cfg.template getParameter<double>("nStubsMin")) {}
167  bool operator()(const L1Track& t) const { return t.getStubRefs().size() >= nStubsMin_; }
168 
169  private:
170  double nStubsMin_;
171  };
175  : nStubsMin_(cfg.template getParameter<double>("nStubsMin")) {}
176  bool operator()(const L1Track& t) const { return t.getNStubs() >= nStubsMin_; }
177 
178  private:
179  double nStubsMin_;
180  };
183  : nPSStubsMin_(nStubsMin), tTopo_(tTopo) {}
185  : nPSStubsMin_(cfg.template getParameter<double>("nPSStubsMin")), tTopo_(tTopo) {}
186  bool operator()(const L1Track& t) const {
187  int nPSStubs = 0;
188  for (const auto& stub : t.getStubRefs()) {
189  DetId detId(stub->getDetId());
190  if (detId.det() == DetId::Detector::Tracker) {
191  if ((detId.subdetId() == StripSubdetector::TOB && tTopo_.tobLayer(detId) <= 3) ||
192  (detId.subdetId() == StripSubdetector::TID && tTopo_.tidRing(detId) <= 9))
193  nPSStubs++;
194  }
195  }
196  return nPSStubs >= nPSStubsMin_;
197  }
198 
199  private:
200  double nPSStubsMin_;
202  };
206  : bendChi2Max_(cfg.template getParameter<double>("reducedBendChi2Max")) {}
207  bool operator()(const L1Track& t) const { return t.stubPtConsistency() < bendChi2Max_; }
208 
209  private:
210  double bendChi2Max_;
211  };
215  : bendChi2Max_(cfg.template getParameter<double>("reducedBendChi2Max")) {}
216  bool operator()(const L1Track& t) const { return t.getBendChi2() < bendChi2Max_; }
217 
218  private:
219  double bendChi2Max_;
220  };
224  : reducedChi2RZMax_(cfg.template getParameter<double>("reducedChi2RZMax")) {}
225  bool operator()(const L1Track& t) const { return t.chi2ZRed() < reducedChi2RZMax_; }
226 
227  private:
229  };
233  : reducedChi2RZMax_(cfg.template getParameter<double>("reducedChi2RZMax")) {}
234  bool operator()(const L1Track& t) const { return t.getChi2RZ() < reducedChi2RZMax_; }
235 
236  private:
238  };
242  : reducedChi2RPhiMax_(cfg.template getParameter<double>("reducedChi2RPhiMax")) {}
243  bool operator()(const L1Track& t) const { return t.chi2XYRed() < reducedChi2RPhiMax_; }
244 
245  private:
247  };
251  : reducedChi2RPhiMax_(cfg.template getParameter<double>("reducedChi2RPhiMax")) {}
252  bool operator()(const L1Track& t) const { return t.getChi2RPhi() < reducedChi2RPhiMax_; }
253 
254  private:
256  };
258  TTTrackDeltaZMaxSelector(const std::vector<double>& deltaZMaxEtaBounds, const std::vector<double>& deltaZMax)
261  : deltaZMaxEtaBounds_(cfg.template getParameter<double>("deltaZMaxEtaBounds")),
262  deltaZMax_(cfg.template getParameter<double>("deltaZMax")) {}
263  bool operator()(const L1Track& t, const l1t::Vertex& v) const {
264  size_t etaIndex =
265  std::lower_bound(deltaZMaxEtaBounds_.begin(), deltaZMaxEtaBounds_.end(), std::abs(t.momentum().eta())) -
266  deltaZMaxEtaBounds_.begin() - 1;
267  if (etaIndex > deltaZMax_.size() - 1)
268  etaIndex = deltaZMax_.size() - 1;
269  return std::abs(v.z0() - t.z0()) <= deltaZMax_[etaIndex];
270  }
271 
272  private:
273  std::vector<double> deltaZMaxEtaBounds_;
274  std::vector<double> deltaZMax_;
275  };
277  TTTrackWordDeltaZMaxSelector(const std::vector<double>& deltaZMaxEtaBounds, const std::vector<double>& deltaZMax)
280  : deltaZMaxEtaBounds_(cfg.template getParameter<double>("deltaZMaxEtaBounds")),
281  deltaZMax_(cfg.template getParameter<double>("deltaZMax")) {}
282  bool operator()(const L1Track& t, const l1t::VertexWord& v) const {
283  size_t etaIndex =
284  std::lower_bound(deltaZMaxEtaBounds_.begin(), deltaZMaxEtaBounds_.end(), std::abs(t.momentum().eta())) -
285  deltaZMaxEtaBounds_.begin() - 1;
286  if (etaIndex > deltaZMax_.size() - 1)
287  etaIndex = deltaZMax_.size() - 1;
288  return std::abs(v.z0() - t.getZ0()) <= deltaZMax_[etaIndex];
289  }
290 
291  private:
292  std::vector<double> deltaZMaxEtaBounds_;
293  std::vector<double> deltaZMax_;
294  };
295 
307 
308  // ----------member data ---------------------------
317  std::vector<double> deltaZMaxEtaBounds_, deltaZMax_;
320  int debug_;
321 };
322 
323 //
324 // constructors and destructor
325 //
327  : l1TracksToken_(consumes<TTTrackCollection>(iConfig.getParameter<edm::InputTag>("l1TracksInputTag"))),
329  outputCollectionName_(iConfig.getParameter<std::string>("outputCollectionName")),
330  cutSet_(iConfig.getParameter<edm::ParameterSet>("cutSet")),
331 
332  ptMin_(cutSet_.getParameter<double>("ptMin")),
333  absEtaMax_(cutSet_.getParameter<double>("absEtaMax")),
334  absZ0Max_(cutSet_.getParameter<double>("absZ0Max")),
335  bendChi2Max_(cutSet_.getParameter<double>("reducedBendChi2Max")),
336  reducedChi2RZMax_(cutSet_.getParameter<double>("reducedChi2RZMax")),
337  reducedChi2RPhiMax_(cutSet_.getParameter<double>("reducedChi2RPhiMax")),
338  nStubsMin_(cutSet_.getParameter<int>("nStubsMin")),
339  nPSStubsMin_(cutSet_.getParameter<int>("nPSStubsMin")),
340  deltaZMaxEtaBounds_(cutSet_.getParameter<std::vector<double>>("deltaZMaxEtaBounds")),
341  deltaZMax_(cutSet_.getParameter<std::vector<double>>("deltaZMax")),
342 
343  useDisplacedTracksDeltaZOverride_(iConfig.getParameter<double>("useDisplacedTracksDeltaZOverride")),
344  processSimulatedTracks_(iConfig.getParameter<bool>("processSimulatedTracks")),
345  processEmulatedTracks_(iConfig.getParameter<bool>("processEmulatedTracks")),
346  debug_(iConfig.getParameter<int>("debug")) {
347  // Confirm the the configuration makes sense
349  throw cms::Exception("You must process at least one of the track collections (simulated or emulated).");
350  }
351 
352  if (deltaZMax_.size() != deltaZMaxEtaBounds_.size() - 1) {
353  throw cms::Exception("The number of deltaZ cuts does not match the number of eta bins!");
354  }
355 
357  deltaZMax_ = std::vector<double>(deltaZMax_.size(), useDisplacedTracksDeltaZOverride_);
358  }
359 
360  // Get additional input tags and define the EDM output based on the previous configuration parameters
361  doDeltaZCutSim_ = false;
362  doDeltaZCutEmu_ = false;
364  produces<TTTrackRefCollection>(outputCollectionName_);
365  if (iConfig.exists("l1VerticesInputTag")) {
366  l1VerticesToken_ = consumes<l1t::VertexCollection>(iConfig.getParameter<edm::InputTag>("l1VerticesInputTag"));
367  doDeltaZCutSim_ = true;
368  produces<TTTrackRefCollection>(outputCollectionName_ + "Associated");
369  }
370  }
372  produces<TTTrackRefCollection>(outputCollectionName_ + "Emulation");
373  if (iConfig.exists("l1VerticesEmulationInputTag")) {
375  consumes<l1t::VertexWordCollection>(iConfig.getParameter<edm::InputTag>("l1VerticesEmulationInputTag"));
376  doDeltaZCutEmu_ = true;
377  produces<TTTrackRefCollection>(outputCollectionName_ + "AssociatedEmulation");
378  }
379  }
380 }
381 
383 
384 //
385 // member functions
386 //
387 
389  const TTTrackRefCollectionUPtr& vTTTrackOutput,
390  const TTTrackRefCollectionUPtr& vTTTrackEmulationOutput,
391  const TTTrackRefCollectionUPtr& vTTTrackAssociatedOutput,
392  const TTTrackRefCollectionUPtr& vTTTrackAssociatedEmulationOutput) const {
393  edm::LogInfo log("L1TrackSelectionProducer");
394  log << "The original track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are ... \n";
395  for (const auto& track : *l1TracksHandle) {
396  printTrackInfo(log, track, debug_ >= 4);
397  }
398  log << "\t---\n\tNumber of tracks in this selection = " << l1TracksHandle->size() << "\n\n";
400  log << "The selected track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are ... \n";
401  for (const auto& track : *vTTTrackOutput) {
402  printTrackInfo(log, *track, debug_ >= 4);
403  }
404  log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackOutput->size() << "\n\n";
405  }
407  log << "The emulation selected track collection (pt, eta, phi, nstub, bendchi2, chi2rz, chi2rphi, z0) values are "
408  "... \n";
409  for (const auto& track : *vTTTrackEmulationOutput) {
410  printTrackInfo(log, *track, debug_ >= 4);
411  }
412  log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackEmulationOutput->size() << "\n\n";
413  }
415  TTTrackRefCollection inSimButNotEmu;
416  TTTrackRefCollection inEmuButNotSim;
417  std::set_difference(vTTTrackOutput->begin(),
418  vTTTrackOutput->end(),
419  vTTTrackEmulationOutput->begin(),
420  vTTTrackEmulationOutput->end(),
421  std::back_inserter(inSimButNotEmu));
422  std::set_difference(vTTTrackEmulationOutput->begin(),
423  vTTTrackEmulationOutput->end(),
424  vTTTrackOutput->begin(),
425  vTTTrackOutput->end(),
426  std::back_inserter(inEmuButNotSim));
427  log << "The set of tracks selected via cuts on the simulated values which are not in the set of tracks selected "
428  "by cutting on the emulated values ... \n";
429  for (const auto& track : inSimButNotEmu) {
430  printTrackInfo(log, *track, debug_ >= 3);
431  }
432  log << "\t---\n\tNumber of tracks in this selection = " << inSimButNotEmu.size() << "\n\n"
433  << "The set of tracks selected via cuts on the emulated values which are not in the set of tracks selected "
434  "by cutting on the simulated values ... \n";
435  for (const auto& track : inEmuButNotSim) {
436  printTrackInfo(log, *track, debug_ >= 3);
437  }
438  log << "\t---\n\tNumber of tracks in this selection = " << inEmuButNotSim.size() << "\n\n";
439  }
441  log << "The selected and leading vertex associated track collection (pt, eta, phi, nstub, bendchi2, chi2rz, "
442  "chi2rphi, z0) values are ... \n";
443  for (const auto& track : *vTTTrackAssociatedOutput) {
444  printTrackInfo(log, *track, debug_ >= 4);
445  }
446  log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackAssociatedOutput->size() << "\n\n";
447  }
449  log << "The emulation selected and leading vertex associated track collection (pt, eta, phi, nstub, bendchi2, "
450  "chi2rz, chi2rphi, z0) values are "
451  "... \n";
452  for (const auto& track : *vTTTrackAssociatedEmulationOutput) {
453  printTrackInfo(log, *track, debug_ >= 4);
454  }
455  log << "\t---\n\tNumber of tracks in this selection = " << vTTTrackAssociatedEmulationOutput->size() << "\n\n";
456  }
458  TTTrackRefCollection inSimButNotEmu;
459  TTTrackRefCollection inEmuButNotSim;
460  std::set_difference(vTTTrackAssociatedOutput->begin(),
461  vTTTrackAssociatedOutput->end(),
462  vTTTrackAssociatedEmulationOutput->begin(),
463  vTTTrackAssociatedEmulationOutput->end(),
464  std::back_inserter(inSimButNotEmu));
465  std::set_difference(vTTTrackAssociatedEmulationOutput->begin(),
466  vTTTrackAssociatedEmulationOutput->end(),
467  vTTTrackAssociatedOutput->begin(),
468  vTTTrackAssociatedOutput->end(),
469  std::back_inserter(inEmuButNotSim));
470  log << "The set of tracks selected via cuts on the simulated values which are not in the set of tracks selected "
471  "by cutting on the emulated values ... \n";
472  for (const auto& track : inSimButNotEmu) {
473  printTrackInfo(log, *track, debug_ >= 3);
474  }
475  log << "\t---\n\tNumber of tracks in this selection = " << inSimButNotEmu.size() << "\n\n"
476  << "The set of tracks selected via cuts on the emulated values which are not in the set of tracks selected "
477  "by cutting on the simulated values ... \n";
478  for (const auto& track : inEmuButNotSim) {
479  printTrackInfo(log, *track, debug_ >= 3);
480  }
481  log << "\t---\n\tNumber of tracks in this selection = " << inEmuButNotSim.size() << "\n\n";
482  }
483 }
484 
485 void L1TrackSelectionProducer::printTrackInfo(edm::LogInfo& log, const L1Track& track, bool printEmulation) const {
486  log << "\t(" << track.momentum().perp() << ", " << track.momentum().eta() << ", " << track.momentum().phi() << ", "
487  << track.getStubRefs().size() << ", " << track.stubPtConsistency() << ", " << track.chi2ZRed() << ", "
488  << track.chi2XYRed() << ", " << track.z0() << ")\n";
489 
490  if (printEmulation) {
491  ap_uint<TrackBitWidths::kPtSize> ptEmulationBits = track.getTrackWord()(
492  TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1, TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
493  ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize> ptEmulation;
494  ptEmulation.V = ptEmulationBits.range();
495  TTTrack_TrackWord::tanl_t etaEmulationBits = track.getTanlWord();
496  ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize> etaEmulation;
497  etaEmulation.V = etaEmulationBits.range();
498  log << "\t\t(" << ptEmulation.to_double() << ", " << etaEmulation.to_double() << ", " << track.getPhi() << ", "
499  << track.getNStubs() << ", " << track.getBendChi2() << ", " << track.getChi2RZ() << ", " << track.getChi2RPhi()
500  << ", " << track.getZ0() << ")\n";
501  }
502 }
503 
504 // ------------ method called to produce the data ------------
506  auto vTTTrackOutput = std::make_unique<TTTrackRefCollection>();
507  auto vTTTrackAssociatedOutput = std::make_unique<TTTrackRefCollection>();
508  auto vTTTrackEmulationOutput = std::make_unique<TTTrackRefCollection>();
509  auto vTTTrackAssociatedEmulationOutput = std::make_unique<TTTrackRefCollection>();
510 
511  // Tracker Topology
512  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
513 
514  TTTrackCollectionHandle l1TracksHandle;
515  edm::Handle<l1t::VertexCollection> l1VerticesHandle;
516  edm::Handle<l1t::VertexWordCollection> l1VerticesEmulationHandle;
517 
518  l1t::Vertex leadingVertex;
519  l1t::VertexWord leadingEmulationVertex;
520 
521  iEvent.getByToken(l1TracksToken_, l1TracksHandle);
522  size_t nOutputApproximate = l1TracksHandle->size();
524  if (doDeltaZCutSim_) {
525  iEvent.getByToken(l1VerticesToken_, l1VerticesHandle);
526  leadingVertex = l1VerticesHandle->at(0);
527  if (debug_ >= 2) {
528  edm::LogInfo("L1TrackSelectionProducer") << "leading vertex z0 = " << leadingVertex.z0();
529  }
530  }
531  vTTTrackOutput->reserve(nOutputApproximate);
532  vTTTrackAssociatedOutput->reserve(nOutputApproximate);
533  }
535  if (doDeltaZCutEmu_) {
536  iEvent.getByToken(l1VerticesEmulationToken_, l1VerticesEmulationHandle);
537  leadingEmulationVertex = l1VerticesEmulationHandle->at(0);
538  if (debug_ >= 2) {
539  edm::LogInfo("L1TrackSelectionProducer") << "leading emulation vertex z0 = " << leadingEmulationVertex.z0();
540  }
541  }
542  vTTTrackEmulationOutput->reserve(nOutputApproximate);
543  vTTTrackAssociatedEmulationOutput->reserve(nOutputApproximate);
544  }
545 
552  TTTrackNPSStubsMinSelector nPSStubsSel(nPSStubsMin_, tTopo);
553 
554  for (size_t i = 0; i < nOutputApproximate; i++) {
555  const auto& track = l1TracksHandle->at(i);
556 
557  // Select tracks based on the floating point TTTrack
558  if (processSimulatedTracks_ && kinSel(track) && nPSStubsSel(track) && chi2Sel(track)) {
559  vTTTrackOutput->push_back(TTTrackRef(l1TracksHandle, i));
560  if (doDeltaZCutSim_ && deltaZSel(track, leadingVertex)) {
561  vTTTrackAssociatedOutput->push_back(TTTrackRef(l1TracksHandle, i));
562  }
563  }
564 
565  // Select tracks based on the bitwise accurate TTTrack_TrackWord
566  if (processEmulatedTracks_ && kinSelEmu(track) && chi2SelEmu(track)) {
567  vTTTrackEmulationOutput->push_back(TTTrackRef(l1TracksHandle, i));
568  if (doDeltaZCutEmu_ && deltaZSelEmu(track, leadingEmulationVertex)) {
569  vTTTrackAssociatedEmulationOutput->push_back(TTTrackRef(l1TracksHandle, i));
570  }
571  }
572  }
573 
574  if (debug_ >= 2) {
575  printDebugInfo(l1TracksHandle,
576  vTTTrackOutput,
577  vTTTrackEmulationOutput,
578  vTTTrackAssociatedOutput,
579  vTTTrackAssociatedEmulationOutput);
580  }
581 
582  // Put the outputs into the event
584  iEvent.put(std::move(vTTTrackOutput), outputCollectionName_);
585  if (doDeltaZCutSim_) {
586  iEvent.put(std::move(vTTTrackAssociatedOutput), outputCollectionName_ + "Associated");
587  }
588  }
590  iEvent.put(std::move(vTTTrackEmulationOutput), outputCollectionName_ + "Emulation");
591  if (doDeltaZCutEmu_) {
592  iEvent.put(std::move(vTTTrackAssociatedEmulationOutput), outputCollectionName_ + "AssociatedEmulation");
593  }
594  }
595 }
596 
597 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
599  //L1TrackSelectionProducer
601  desc.add<edm::InputTag>("l1TracksInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
602  desc.addOptional<edm::InputTag>("l1VerticesInputTag", edm::InputTag("L1VertexFinder", "l1vertices"));
603  desc.addOptional<edm::InputTag>("l1VerticesEmulationInputTag",
604  edm::InputTag("L1VertexFinderEmulator", "l1verticesEmulation"));
605  desc.add<std::string>("outputCollectionName", "Level1TTTracksSelected");
606  {
607  edm::ParameterSetDescription descCutSet;
608  descCutSet.add<double>("ptMin", 2.0)->setComment("pt must be greater than this value, [GeV]");
609  descCutSet.add<double>("absEtaMax", 2.4)->setComment("absolute value of eta must be less than this value");
610  descCutSet.add<double>("absZ0Max", 15.0)->setComment("z0 must be less than this value, [cm]");
611  descCutSet.add<int>("nStubsMin", 4)->setComment("number of stubs must be greater than or equal to this value");
612  descCutSet.add<int>("nPSStubsMin", 0)
613  ->setComment("number of stubs in the PS Modules must be greater than or equal to this value");
614 
615  descCutSet.add<double>("reducedBendChi2Max", 2.25)->setComment("bend chi2 must be less than this value");
616  descCutSet.add<double>("reducedChi2RZMax", 5.0)->setComment("chi2rz/dof must be less than this value");
617  descCutSet.add<double>("reducedChi2RPhiMax", 20.0)->setComment("chi2rphi/dof must be less than this value");
618 
619  descCutSet.add<std::vector<double>>("deltaZMaxEtaBounds", {0.0, 0.7, 1.0, 1.2, 1.6, 2.0, 2.4})
620  ->setComment("these values define the bin boundaries in |eta|");
621  descCutSet.add<std::vector<double>>("deltaZMax", {0.37, 0.50, 0.60, 0.75, 1.00, 1.60})
622  ->setComment(
623  "delta z must be less than these values, there will be one less value here than in deltaZMaxEtaBounds, "
624  "[cm]");
625  desc.add<edm::ParameterSetDescription>("cutSet", descCutSet);
626  }
627  desc.add<double>("useDisplacedTracksDeltaZOverride", -1.0)
628  ->setComment("override the deltaZ cut value for displaced tracks");
629  desc.add<bool>("processSimulatedTracks", true)
630  ->setComment("return selected tracks after cutting on the floating point values");
631  desc.add<bool>("processEmulatedTracks", true)
632  ->setComment("return selected tracks after cutting on the bitwise emulated values");
633  desc.add<int>("debug", 0)->setComment("Verbosity levels: 0, 1, 2, 3");
634  descriptions.addWithDefaultLabel(desc);
635 }
636 
637 //define this as a plug-in
L1TrackSelectionProducer(const edm::ParameterSet &)
void setComment(std::string const &value)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
unsigned int tobLayer(const DetId &id) const
edm::Handle< TTTrackCollection > TTTrackCollectionHandle
TTTrackDeltaZMaxSelector(const std::vector< double > &deltaZMaxEtaBounds, const std::vector< double > &deltaZMax)
AndSelector< TTTrackBendChi2MaxSelector, TTTrackChi2RZMaxSelector, TTTrackChi2RPhiMaxSelector > TTTrackBendChi2Chi2RZChi2RPhiMaxSelector
AndSelector< TTTrackWordBendChi2MaxSelector, TTTrackWordChi2RZMaxSelector, TTTrackWordChi2RPhiMaxSelector > TTTrackWordBendChi2Chi2RZChi2RPhiMaxSelector
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
AndSelector< TTTrackWordPtMinSelector, TTTrackWordAbsEtaMaxSelector, TTTrackWordAbsZ0MaxSelector, TTTrackWordNStubsMinSelector > TTTrackWordPtMinEtaMaxZ0MaxNStubsMinSelector
void printDebugInfo(const TTTrackCollectionHandle &l1TracksHandle, const TTTrackRefCollectionUPtr &vTTTrackOutput, const TTTrackRefCollectionUPtr &vTTTrackEmulationOutput, const TTTrackRefCollectionUPtr &vTTTrackAssociatedOutput, const TTTrackRefCollectionUPtr &vTTTrackAssociatedEmulationOutput) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
AndSelector< TTTrackPtMinSelector, TTTrackAbsEtaMaxSelector, TTTrackAbsZ0MaxSelector, TTTrackNStubsMinSelector > TTTrackPtMinEtaMaxZ0MaxNStubsMinSelector
constexpr float ptMin
TTTrackNPSStubsMinSelector(const edm::ParameterSet &cfg, const TrackerTopology &tTopo)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::unique_ptr< TTTrackRefCollection > TTTrackRefCollectionUPtr
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
TTTrackNPSStubsMinSelector(double nStubsMin, const TrackerTopology &tTopo)
int iEvent
Definition: GenABIO.cc:224
std::vector< double > deltaZMaxEtaBounds_
edm::EDGetTokenT< l1t::VertexCollection > l1VerticesToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Ref< TTTrackCollection > TTTrackRef
std::vector< TTTrackRef > TTTrackRefCollection
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
bool getData(T &iHolder) const
Definition: EventSetup.h:122
bool operator()(const L1Track &t, const l1t::Vertex &v) const
const edm::ParameterSet cutSet_
static constexpr auto TOB
TTTrack< Ref_Phase2TrackerDigi_ > L1Track
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< L1Track > TTTrackCollection
Log< level::Info, false > LogInfo
Definition: DetId.h:17
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
void printTrackInfo(edm::LogInfo &log, const L1Track &track, bool printEmulation=false) const
double z0() const
Definition: VertexWord.h:156
bool operator()(const L1Track &t, const l1t::VertexWord &v) const
edm::EDGetTokenT< l1t::VertexWordCollection > l1VerticesEmulationToken_
HLT enums.
ap_uint< TrackBitWidths::kTanlSize > tanl_t
unsigned int tidRing(const DetId &id) const
float z0() const
Definition: Vertex.h:21
TTTrackWordDeltaZMaxSelector(const std::vector< double > &deltaZMaxEtaBounds, const std::vector< double > &deltaZMax)
static constexpr auto TID
std::vector< std::string > set_difference(std::vector< std::string > const &v1, std::vector< std::string > const &v2)
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< TTTrackCollection > l1TracksToken_