CMS 3D CMS Logo

DTTrigPhase2Prod.cc
Go to the documentation of this file.
10 
17 
21 
24 
47 
61 
62 // DT trigger GeomUtils
64 
65 //RPC TP
70 
71 #include <fstream>
72 #include <iostream>
73 #include <queue>
74 #include <cmath>
75 
76 using namespace edm;
77 using namespace std;
78 using namespace cmsdt;
79 
81  typedef std::map<DTChamberId, DTDigiCollection, std::less<DTChamberId>> DTDigiMap;
82  typedef DTDigiMap::iterator DTDigiMap_iterator;
83  typedef DTDigiMap::const_iterator DTDigiMap_const_iterator;
84 
85 public:
88 
90  ~DTTrigPhase2Prod() override;
91 
93  void beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
94 
96  void produce(edm::Event& iEvent, const edm::EventSetup& iEventSetup) override;
97 
99  void endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
100 
101  // Methods
102  int rango(const metaPrimitive& mp) const;
103  bool outer(const metaPrimitive& mp) const;
104  bool inner(const metaPrimitive& mp) const;
105  void printmP(const std::string& ss, const metaPrimitive& mP) const;
106  void printmP(const metaPrimitive& mP) const;
107  void printmPC(const std::string& ss, const metaPrimitive& mP) const;
108  void printmPC(const metaPrimitive& mP) const;
109  bool hasPosRF(int wh, int sec) const;
110 
111  // Getter-methods
112  MP_QUALITY getMinimumQuality(void);
113 
114  // Setter-methods
115  void setChiSquareThreshold(float ch2Thr);
116  void setMinimumQuality(MP_QUALITY q);
117 
118  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
119 
120  // data-members
123  std::vector<std::pair<int, MuonPath>> primitives_;
124 
125 private:
126  // Trigger Configuration Manager CCB validity flag
128 
129  // BX offset used to correct DTTPG output
131 
132  // Debug Flag
133  bool debug_;
134  bool dump_;
139 
148 
149  // ParameterSet
152 
153  // Grouping attributes and methods
154  int algo_; // Grouping code
155  std::unique_ptr<MotherGrouping> grouping_obj_;
156  std::unique_ptr<MuonPathAnalyzer> mpathanalyzer_;
157  std::unique_ptr<LateralityProvider> latprovider_;
158  std::unique_ptr<MPFilter> mpathqualityenhancer_;
159  std::unique_ptr<MPFilter> mpathqualityenhancerbayes_;
160  std::unique_ptr<MPFilter> mpathredundantfilter_;
161  std::unique_ptr<MPFilter> mpathhitsfilter_;
162  std::unique_ptr<MuonPathAnalyzer> mpathassociator_;
163  std::unique_ptr<MuonPathConfirmator> mpathconfirmator_;
164  std::unique_ptr<MPFilter> mpathcorfilter_;
165  std::shared_ptr<GlobalCoordsObtainer> globalcoordsobtainer_;
166 
167  // Buffering
171  std::vector<DTDigiCollection*> distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ);
172  void processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
173  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec);
174 
175  // RPC
176  std::unique_ptr<RPCIntegrator> rpc_integrator_;
177  bool useRPC_;
178 
179  void assignIndex(std::vector<metaPrimitive>& inMPaths);
180  void assignIndexPerBX(std::vector<metaPrimitive>& inMPaths);
181  int assignQualityOrder(const metaPrimitive& mP) const;
182 
183  const std::unordered_map<int, int> qmap_;
184 };
185 
186 namespace {
187  struct {
188  bool operator()(std::pair<DTLayerId, DTDigi> a, std::pair<DTLayerId, DTDigi> b) const {
189  return (a.second.time() < b.second.time());
190  }
191  } const DigiTimeOrdering;
192 } // namespace
193 
195  : qmap_({{8, 8}, {7, 7}, {6, 6}, {4, 4}, {3, 3}, {2, 2}, {1, 1}}) {
196  produces<L1Phase2MuDTPhContainer>();
197  produces<L1Phase2MuDTThContainer>();
198  produces<L1Phase2MuDTExtPhContainer>();
199  produces<L1Phase2MuDTExtThContainer>();
200 
201  debug_ = pset.getUntrackedParameter<bool>("debug");
202  dump_ = pset.getUntrackedParameter<bool>("dump");
203 
204  scenario_ = pset.getParameter<int>("scenario");
205 
206  df_extended_ = pset.getParameter<int>("df_extended");
207  max_index_ = pset.getParameter<int>("max_primitives") - 1;
208 
209  dtDigisToken_ = consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiTag"));
210 
211  rpcRecHitsLabel_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("rpcRecHits"));
212  useRPC_ = pset.getParameter<bool>("useRPC");
213 
214  // Choosing grouping scheme:
215  algo_ = pset.getParameter<int>("algo");
216 
217  // shortcuts
218 
219  output_mixer_ = pset.getParameter<bool>("output_mixer");
220  output_latpredictor_ = pset.getParameter<bool>("output_latpredictor");
221  output_slfitter_ = pset.getParameter<bool>("output_slfitter");
222  output_slfilter_ = pset.getParameter<bool>("output_slfilter");
223  output_confirmed_ = pset.getParameter<bool>("output_confirmed");
224  output_matcher_ = pset.getParameter<bool>("output_matcher");
225  allow_confirmation_ = pset.getParameter<bool>("allow_confirmation");
226 
227  edm::ConsumesCollector consumesColl(consumesCollector());
228  globalcoordsobtainer_ = std::make_shared<GlobalCoordsObtainer>(pset);
229  globalcoordsobtainer_->generate_luts();
230 
231  if (algo_ == PseudoBayes) {
232  grouping_obj_ =
233  std::make_unique<PseudoBayesGrouping>(pset.getParameter<edm::ParameterSet>("PseudoBayesPattern"), consumesColl);
234  } else if (algo_ == HoughTrans) {
235  grouping_obj_ =
236  std::make_unique<HoughGrouping>(pset.getParameter<edm::ParameterSet>("HoughGrouping"), consumesColl);
237  } else {
238  grouping_obj_ = std::make_unique<TrapezoidalGrouping>(pset, consumesColl);
239  }
240 
241  if (algo_ == Standard) {
242  if (debug_)
243  LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: JM analyzer";
244  mpathanalyzer_ = std::make_unique<MuonPathSLFitter>(pset, consumesColl, globalcoordsobtainer_);
245  latprovider_ = std::make_unique<LateralityCoarsedProvider>(pset, consumesColl);
246  } else {
247  if (debug_)
248  LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: Full chamber analyzer";
249  mpathanalyzer_ = std::make_unique<MuonPathAnalyzerInChamber>(pset, consumesColl, globalcoordsobtainer_);
250  }
251 
252  // Getting buffer option
253  activateBuffer_ = pset.getParameter<bool>("activateBuffer");
254  superCellhalfspacewidth_ = pset.getParameter<int>("superCellspacewidth") / 2;
255  superCelltimewidth_ = pset.getParameter<double>("superCelltimewidth");
256 
257  mpathqualityenhancer_ = std::make_unique<MPSLFilter>(pset);
258  mpathqualityenhancerbayes_ = std::make_unique<MPQualityEnhancerFilterBayes>(pset);
259  mpathredundantfilter_ = std::make_unique<MPRedundantFilter>(pset);
260  mpathhitsfilter_ = std::make_unique<MPCleanHitsFilter>(pset);
261  mpathconfirmator_ = std::make_unique<MuonPathConfirmator>(pset, consumesColl);
262  mpathassociator_ = std::make_unique<MuonPathCorFitter>(pset, consumesColl, globalcoordsobtainer_);
263  mpathcorfilter_ = std::make_unique<MPCorFilter>(pset);
264  rpc_integrator_ = std::make_unique<RPCIntegrator>(pset, consumesColl);
265 
266  dtGeomH = esConsumes<DTGeometry, MuonGeometryRecord, edm::Transition::BeginRun>();
267 }
268 
270  if (debug_)
271  LogDebug("DTTrigPhase2Prod") << "DTp2: calling destructor" << std::endl;
272 }
273 
274 void DTTrigPhase2Prod::beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
275  if (debug_)
276  LogDebug("DTTrigPhase2Prod") << "beginRun " << iRun.id().run();
277  if (debug_)
278  LogDebug("DTTrigPhase2Prod") << "beginRun: getting DT geometry";
279 
280  grouping_obj_->initialise(iEventSetup); // Grouping object initialisation
281  mpathanalyzer_->initialise(iEventSetup); // Analyzer object initialisation
282  mpathqualityenhancer_->initialise(iEventSetup); // Filter object initialisation
283  mpathredundantfilter_->initialise(iEventSetup); // Filter object initialisation
284  mpathqualityenhancerbayes_->initialise(iEventSetup); // Filter object initialisation
285  mpathhitsfilter_->initialise(iEventSetup);
286  mpathassociator_->initialise(iEventSetup); // Associator object initialisation
287  mpathcorfilter_->initialise(iEventSetup);
288 
289  if (auto geom = iEventSetup.getHandle(dtGeomH)) {
290  dtGeo_ = &(*geom);
291  }
292 }
293 
294 void DTTrigPhase2Prod::produce(Event& iEvent, const EventSetup& iEventSetup) {
295  if (debug_)
296  LogDebug("DTTrigPhase2Prod") << "produce";
298  iEvent.getByToken(dtDigisToken_, dtdigis);
299 
300  if (debug_)
301  LogDebug("DTTrigPhase2Prod") << "\t Getting the RPC RecHits" << std::endl;
303  iEvent.getByToken(rpcRecHitsLabel_, rpcRecHits);
304 
306  // GROUPING CODE:
308 
309  DTDigiMap digiMap;
311  for (const auto& detUnitIt : *dtdigis) {
312  const DTLayerId& layId = detUnitIt.first;
313  const DTChamberId chambId = layId.superlayerId().chamberId();
314  const DTDigiCollection::Range& range = detUnitIt.second;
315  digiMap[chambId].put(range, layId);
316  }
317 
318  // generate a list muon paths for each event!!!
319  if (debug_ && activateBuffer_)
320  LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber using a buffer and super cells.";
321  else if (debug_)
322  LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber.";
323 
324  std::map<int, MuonPathPtrs> muonpaths;
325  for (const auto& ich : dtGeo_->chambers()) {
326  // The code inside this for loop would ideally later fit inside a trigger unit (in principle, a DT station) of the future Phase 2 DT Trigger.
327  const DTChamber* chamb = ich;
328  DTChamberId chid = chamb->id();
329  DTDigiMap_iterator dmit = digiMap.find(chid);
330 
331  if (dmit == digiMap.end())
332  continue;
333 
334  if (activateBuffer_) { // Use buffering (per chamber) or not
335  // Import digis from the station
336  std::vector<std::pair<DTLayerId, DTDigi>> tmpvec;
337  tmpvec.clear();
338 
339  for (const auto& dtLayerIdIt : (*dmit).second) {
340  for (DTDigiCollection::const_iterator digiIt = (dtLayerIdIt.second).first;
341  digiIt != (dtLayerIdIt.second).second;
342  digiIt++) {
343  tmpvec.emplace_back(dtLayerIdIt.first, *digiIt);
344  }
345  }
346 
347  // Check to enhance CPU time usage
348  if (tmpvec.empty())
349  continue;
350 
351  // Order digis depending on TDC time and insert them into a queue (FIFO buffer). TODO: adapt for MC simulations.
352  std::sort(tmpvec.begin(), tmpvec.end(), DigiTimeOrdering);
353  std::queue<std::pair<DTLayerId, DTDigi>> timequeue;
354 
355  for (const auto& elem : tmpvec)
356  timequeue.emplace(elem);
357  tmpvec.clear();
358 
359  // Distribute the digis from the queue into supercells
360  std::vector<DTDigiCollection*> superCells;
361  superCells = distribDigis(timequeue);
362 
363  // Process each supercell & collect the resulting muonpaths (as the muonpaths std::vector is only enlarged each time
364  // the groupings access it, it's not needed to "collect" the final products).
365 
366  while (!superCells.empty()) {
367  grouping_obj_->run(iEvent, iEventSetup, *(superCells.back()), muonpaths[chid.rawId()]);
368  superCells.pop_back();
369  }
370  } else {
371  grouping_obj_->run(iEvent, iEventSetup, (*dmit).second, muonpaths[chid.rawId()]);
372  }
373  }
374  digiMap.clear();
375 
376  if (dump_) {
377  for (auto& ch_muonpaths : muonpaths) {
378  for (unsigned int i = 0; i < ch_muonpaths.second.size(); i++) {
379  stringstream ss;
380  ss << iEvent.id().event() << " mpath " << i << ": ";
381  for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
382  ss << ch_muonpaths.second.at(i)->primitive(lay)->channelId() << " ";
383  for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
384  ss << ch_muonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
385  for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
386  ss << ch_muonpaths.second.at(i)->primitive(lay)->laterality() << " ";
387  LogInfo("DTTrigPhase2Prod") << ss.str();
388  }
389  }
390  }
391 
392  std::map<int, std::vector<lat_vector>> lateralities;
393  if (!output_mixer_) {
394  for (auto& ch_muonpaths : muonpaths) {
395  if (algo_ == Standard) {
396  latprovider_->run(iEvent, iEventSetup, ch_muonpaths.second, lateralities[ch_muonpaths.first]);
397  }
398  }
399  }
400 
401  // FILTER GROUPING
402  std::map<int, MuonPathPtrs> filteredmuonpaths;
403  for (auto& ch_muonpaths : muonpaths) {
404  if (algo_ == Standard) {
405  mpathredundantfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
406  } else {
407  mpathhitsfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
408  }
409  }
410 
411  if (dump_) {
412  for (auto& ch_filteredmuonpaths : filteredmuonpaths) {
413  for (unsigned int i = 0; i < ch_filteredmuonpaths.second.size(); i++) {
414  stringstream ss;
415  ss << iEvent.id().event() << " filt. mpath " << i << ": ";
416  for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
417  ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->channelId() << " ";
418  for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
419  ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
420  LogInfo("DTTrigPhase2Prod") << ss.str();
421  }
422  }
423  }
424 
426 
430 
431  if (debug_) {
432  for (auto& ch_muonpaths : muonpaths) {
433  LogDebug("DTTrigPhase2Prod") << "MUON PATHS found: " << ch_muonpaths.second.size() << " ("
434  << filteredmuonpaths[ch_muonpaths.first].size() << ") in event "
435  << iEvent.id().event();
436  }
437  }
438  if (debug_)
439  LogDebug("DTTrigPhase2Prod") << "filling NmetaPrimtives" << std::endl;
440  std::map<int, std::vector<metaPrimitive>> metaPrimitives;
441  std::map<int, MuonPathPtrs> outmpaths;
442  if (algo_ == Standard) {
443  if (debug_)
444  LogDebug("DTTrigPhase2Prod") << "Fitting 1SL ";
445  for (auto& ch_muonpaths : muonpaths) { // FIXME, do we need filtered muonpaths?
447  mpathanalyzer_->run(iEvent,
448  iEventSetup,
449  ch_muonpaths.second,
450  lateralities[ch_muonpaths.first],
451  metaPrimitives[ch_muonpaths.first]);
452  else if (output_mixer_) {
453  for (auto& inMPath : ch_muonpaths.second) {
454  auto sl = inMPath->primitive(0)->superLayerId(); // 0, 1, 2
455  int selected_lay = 1;
456  if (inMPath->primitive(0)->tdcTimeStamp() != -1)
457  selected_lay = 0;
458  int dumLayId = inMPath->primitive(selected_lay)->cameraId();
459  auto dtDumlayerId = DTLayerId(dumLayId);
460  DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
461  if (sl == 0)
462  metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
463  -1,
464  -1,
465  -1,
466  -1,
467  -1,
468  -1,
469  -1,
470  -1,
471  -1,
472  inMPath->primitive(0)->channelId(),
473  inMPath->primitive(0)->tdcTimeStamp(),
474  -1,
475  inMPath->primitive(1)->channelId(),
476  inMPath->primitive(1)->tdcTimeStamp(),
477  -1,
478  inMPath->primitive(2)->channelId(),
479  inMPath->primitive(2)->tdcTimeStamp(),
480  -1,
481  inMPath->primitive(3)->channelId(),
482  inMPath->primitive(3)->tdcTimeStamp(),
483  -1,
484  -1,
485  -1,
486  -1,
487  -1,
488  -1,
489  -1,
490  -1,
491  -1,
492  -1,
493  -1,
494  -1,
495  -1,
496  -1}));
497  else
498  metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
499  -1,
500  -1,
501  -1,
502  -1,
503  -1,
504  -1,
505  -1,
506  -1,
507  -1,
508  -1,
509  -1,
510  -1,
511  -1,
512  -1,
513  -1,
514  -1,
515  -1,
516  -1,
517  -1,
518  -1,
519  -1,
520  inMPath->primitive(0)->channelId(),
521  inMPath->primitive(0)->tdcTimeStamp(),
522  -1,
523  inMPath->primitive(1)->channelId(),
524  inMPath->primitive(1)->tdcTimeStamp(),
525  -1,
526  inMPath->primitive(2)->channelId(),
527  inMPath->primitive(2)->tdcTimeStamp(),
528  -1,
529  inMPath->primitive(3)->channelId(),
530  inMPath->primitive(3)->tdcTimeStamp(),
531  -1,
532  -1}));
533  }
534  } else if (output_latpredictor_) {
535  int imp = -1;
536  for (auto& inMPath : ch_muonpaths.second) {
537  imp++;
538  auto sl = inMPath->primitive(0)->superLayerId(); // 0, 1, 2
539  int selected_lay = 1;
540  if (inMPath->primitive(0)->tdcTimeStamp() != -1)
541  selected_lay = 0;
542  int dumLayId = inMPath->primitive(selected_lay)->cameraId();
543  auto dtDumlayerId = DTLayerId(dumLayId);
544  DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
545  for (auto& latcomb : lateralities[ch_muonpaths.first][imp]) {
546  if (sl == 0)
547  metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
548  -1,
549  -1,
550  -1,
551  -1,
552  -1,
553  -1,
554  -1,
555  -1,
556  -1,
557  inMPath->primitive(0)->channelId(),
558  inMPath->primitive(0)->tdcTimeStamp(),
559  latcomb[0],
560  inMPath->primitive(1)->channelId(),
561  inMPath->primitive(1)->tdcTimeStamp(),
562  latcomb[1],
563  inMPath->primitive(2)->channelId(),
564  inMPath->primitive(2)->tdcTimeStamp(),
565  latcomb[2],
566  inMPath->primitive(3)->channelId(),
567  inMPath->primitive(3)->tdcTimeStamp(),
568  latcomb[3],
569  -1,
570  -1,
571  -1,
572  -1,
573  -1,
574  -1,
575  -1,
576  -1,
577  -1,
578  -1,
579  -1,
580  -1,
581  -1}));
582  else
583  metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
584  -1,
585  -1,
586  -1,
587  -1,
588  -1,
589  -1,
590  -1,
591  -1,
592  -1,
593  -1,
594  -1,
595  -1,
596  -1,
597  -1,
598  -1,
599  -1,
600  -1,
601  -1,
602  -1,
603  -1,
604  -1,
605  inMPath->primitive(0)->channelId(),
606  inMPath->primitive(0)->tdcTimeStamp(),
607  latcomb[0],
608  inMPath->primitive(1)->channelId(),
609  inMPath->primitive(1)->tdcTimeStamp(),
610  latcomb[1],
611  inMPath->primitive(2)->channelId(),
612  inMPath->primitive(2)->tdcTimeStamp(),
613  latcomb[2],
614  inMPath->primitive(3)->channelId(),
615  inMPath->primitive(3)->tdcTimeStamp(),
616  latcomb[3],
617  -1}));
618  }
619  }
620  }
621  }
622  } else {
623  // implementation for advanced (2SL) grouping, no filter required..
624  if (debug_)
625  LogDebug("DTTrigPhase2Prod") << "Fitting 2SL at once ";
626  for (auto& ch_muonpaths : muonpaths) {
627  mpathanalyzer_->run(iEvent, iEventSetup, ch_muonpaths.second, outmpaths[ch_muonpaths.first]);
628  }
629  }
630 
632 
633  if (dump_) {
634  for (auto& ch_outmpaths : outmpaths) {
635  for (unsigned int i = 0; i < ch_outmpaths.second.size(); i++) {
636  LogInfo("DTTrigPhase2Prod") << iEvent.id().event() << " mp " << i << ": "
637  << ch_outmpaths.second.at(i)->bxTimeValue() << " "
638  << ch_outmpaths.second.at(i)->horizPos() << " "
639  << ch_outmpaths.second.at(i)->tanPhi() << " " << ch_outmpaths.second.at(i)->phi()
640  << " " << ch_outmpaths.second.at(i)->phiB() << " "
641  << ch_outmpaths.second.at(i)->quality() << " "
642  << ch_outmpaths.second.at(i)->chiSquare();
643  }
644  }
645  for (auto& ch_metaPrimitives : metaPrimitives) {
646  for (unsigned int i = 0; i < ch_metaPrimitives.second.size(); i++) {
647  stringstream ss;
648  ss << iEvent.id().event() << " mp " << i << ": ";
649  printmP(ss.str(), ch_metaPrimitives.second.at(i));
650  }
651  }
652  }
653 
654  muonpaths.clear();
655  filteredmuonpaths.clear();
656 
660 
661  std::map<int, std::vector<metaPrimitive>> confirmedMetaPrimitives;
662  for (auto& ch_metaPrimitives : metaPrimitives) {
664  mpathconfirmator_->run(
665  iEvent, iEventSetup, ch_metaPrimitives.second, dtdigis, confirmedMetaPrimitives[ch_metaPrimitives.first]);
666  else
667  for (auto& mp : ch_metaPrimitives.second) {
668  confirmedMetaPrimitives[ch_metaPrimitives.first].push_back(mp);
669  }
670  }
671 
672  metaPrimitives.clear();
674 
676  // FILTER SECTIONS:
678 
679  if (debug_)
680  LogDebug("DTTrigPhase2Prod") << "declaring new vector for filtered" << std::endl;
681 
682  std::map<int, std::vector<metaPrimitive>> filteredMetaPrimitives;
683  if (algo_ == Standard)
684  for (auto& ch_confirmedMetaPrimitives : confirmedMetaPrimitives) {
685  if (!skip_processing_)
687  iEventSetup,
688  ch_confirmedMetaPrimitives.second,
689  filteredMetaPrimitives[ch_confirmedMetaPrimitives.first]);
690  else
691  for (auto& mp : ch_confirmedMetaPrimitives.second) {
692  filteredMetaPrimitives[ch_confirmedMetaPrimitives.first].push_back(mp);
693  }
694  }
695  if (dump_) {
696  for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
697  for (unsigned int i = 0; i < ch_filteredMetaPrimitives.second.size(); i++) {
698  stringstream ss;
699  ss << iEvent.id().event() << " filtered mp " << i << ": ";
700  printmP(ss.str(), ch_filteredMetaPrimitives.second.at(i));
701  }
702  }
703  }
704 
706  confirmedMetaPrimitives.clear();
707 
708  if (debug_)
709  for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
710  LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
711  << ch_filteredMetaPrimitives.second.size() << " filteredMetaPrimitives (superlayer)"
712  << std::endl;
713  }
714  if (debug_)
715  LogDebug("DTTrigPhase2Prod") << "filteredMetaPrimitives: starting correlations" << std::endl;
716 
720 
721  std::map<int, std::vector<metaPrimitive>> correlatedMetaPrimitives;
722  if (algo_ == Standard) {
723  for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
724  if (!skip_processing_)
726  iEventSetup,
727  ch_filteredMetaPrimitives.second,
728  correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
729  else
730  for (auto& mp : ch_filteredMetaPrimitives.second) {
731  correlatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
732  }
733  }
734  } else {
735  for (auto& ch_outmpaths : outmpaths) {
736  for (const auto& muonpath : ch_outmpaths.second) {
737  correlatedMetaPrimitives[ch_outmpaths.first].emplace_back(muonpath->rawId(),
738  (double)muonpath->bxTimeValue(),
739  muonpath->horizPos(),
740  muonpath->tanPhi(),
741  muonpath->phi(),
742  muonpath->phiB(),
743  muonpath->phi_cmssw(),
744  muonpath->phiB_cmssw(),
745  muonpath->chiSquare(),
746  (int)muonpath->quality(),
747  muonpath->primitive(0)->channelId(),
748  muonpath->primitive(0)->tdcTimeStamp(),
749  muonpath->primitive(0)->laterality(),
750  muonpath->primitive(1)->channelId(),
751  muonpath->primitive(1)->tdcTimeStamp(),
752  muonpath->primitive(1)->laterality(),
753  muonpath->primitive(2)->channelId(),
754  muonpath->primitive(2)->tdcTimeStamp(),
755  muonpath->primitive(2)->laterality(),
756  muonpath->primitive(3)->channelId(),
757  muonpath->primitive(3)->tdcTimeStamp(),
758  muonpath->primitive(3)->laterality(),
759  muonpath->primitive(4)->channelId(),
760  muonpath->primitive(4)->tdcTimeStamp(),
761  muonpath->primitive(4)->laterality(),
762  muonpath->primitive(5)->channelId(),
763  muonpath->primitive(5)->tdcTimeStamp(),
764  muonpath->primitive(5)->laterality(),
765  muonpath->primitive(6)->channelId(),
766  muonpath->primitive(6)->tdcTimeStamp(),
767  muonpath->primitive(6)->laterality(),
768  muonpath->primitive(7)->channelId(),
769  muonpath->primitive(7)->tdcTimeStamp(),
770  muonpath->primitive(7)->laterality());
771  }
772  }
773  }
774 
776 
777  if (debug_)
778  for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
779  LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
780  << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
781  }
782  if (dump_) {
783  for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
784  LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
785  << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
786  }
787  for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
788  for (unsigned int i = 0; i < ch_correlatedMetaPrimitives.second.size(); i++) {
789  stringstream ss;
790  ss << iEvent.id().event() << " correlated mp " << i << ": ";
791  printmPC(ss.str(), ch_correlatedMetaPrimitives.second.at(i));
792  }
793  }
794  }
795 
796  // Correlated Filtering
797  std::map<int, std::vector<metaPrimitive>> filtCorrelatedMetaPrimitives;
798  if (algo_ == Standard) {
799  for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
800  if (!skip_processing_)
801  mpathcorfilter_->run(iEvent,
802  iEventSetup,
803  ch_filteredMetaPrimitives.second,
804  correlatedMetaPrimitives[ch_filteredMetaPrimitives.first],
805  filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
806  else {
807  for (auto& mp : ch_filteredMetaPrimitives.second) {
808  filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
809  }
810  if (output_matcher_)
811  for (auto& mp : correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]) {
812  filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
813  }
814  }
815  }
816  }
817 
818  correlatedMetaPrimitives.clear();
819  filteredMetaPrimitives.clear();
820 
821  double shift_back = 0;
822  if (scenario_ == MC) //scope for MC
823  shift_back = 400;
824  else if (scenario_ == DATA) //scope for data
825  shift_back = 0;
826  else if (scenario_ == SLICE_TEST) //scope for slice test
827  shift_back = 400;
828 
829  // RPC integration
830  if (useRPC_) {
831  rpc_integrator_->initialise(iEventSetup, shift_back);
832  rpc_integrator_->prepareMetaPrimitives(rpcRecHits);
833  for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
834  rpc_integrator_->matchWithDTAndUseRPCTime(ch_correlatedMetaPrimitives.second); // Probably this is a FIXME
835  }
836  rpc_integrator_->makeRPCOnlySegments();
837  rpc_integrator_->storeRPCSingleHits();
838  rpc_integrator_->removeRPCHitsUsed();
839  }
840 
842  vector<L1Phase2MuDTPhDigi> outP2Ph;
843  vector<L1Phase2MuDTExtPhDigi> outExtP2Ph;
844  vector<L1Phase2MuDTThDigi> outP2Th;
845  vector<L1Phase2MuDTExtThDigi> outExtP2Th;
846 
847  // Assigning index value
848  if (!skip_processing_)
849  for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
850  assignIndex(ch_correlatedMetaPrimitives.second);
851  }
852 
853  for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
854  for (const auto& metaPrimitiveIt : ch_correlatedMetaPrimitives.second) {
855  DTChamberId chId(metaPrimitiveIt.rawId);
856  DTSuperLayerId slId(metaPrimitiveIt.rawId);
857  if (debug_)
858  LogDebug("DTTrigPhase2Prod") << "looping in final vector: SuperLayerId" << chId << " x=" << metaPrimitiveIt.x
859  << " quality=" << metaPrimitiveIt.quality
860  << " BX=" << round(metaPrimitiveIt.t0 / 25.) << " index=" << metaPrimitiveIt.index;
861 
862  int sectorTP = chId.sector();
863  //sectors 13 and 14 exist only for the outermost stations for sectors 4 and 10 respectively
864  //due to the larger MB4 that are divided into two.
865  if (sectorTP == 13)
866  sectorTP = 4;
867  if (sectorTP == 14)
868  sectorTP = 10;
869  sectorTP = sectorTP - 1;
870  int sl = 0;
871  if (metaPrimitiveIt.quality < LOWLOWQ || metaPrimitiveIt.quality == CHIGHQ) {
872  if (inner(metaPrimitiveIt))
873  sl = 1;
874  else
875  sl = 3;
876  }
877 
878  float tp_t0 =
879  (metaPrimitiveIt.t0 - shift_back * LHC_CLK_FREQ) * ((float)TIME_TO_TDC_COUNTS / (float)LHC_CLK_FREQ);
880 
881  if (debug_)
882  LogDebug("DTTrigPhase2Prod") << "pushing back phase-2 dataformat carlo-federica dataformat";
883 
884  if (slId.superLayer() != 2) {
885  if (df_extended_ == 1 || df_extended_ == 2) {
886  int pathWireId[8] = {metaPrimitiveIt.wi1,
887  metaPrimitiveIt.wi2,
888  metaPrimitiveIt.wi3,
889  metaPrimitiveIt.wi4,
890  metaPrimitiveIt.wi5,
891  metaPrimitiveIt.wi6,
892  metaPrimitiveIt.wi7,
893  metaPrimitiveIt.wi8};
894 
895  int pathTDC[8] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
896  max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
897  max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
898  max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1),
899  max((int)round(metaPrimitiveIt.tdc5 - shift_back * LHC_CLK_FREQ), -1),
900  max((int)round(metaPrimitiveIt.tdc6 - shift_back * LHC_CLK_FREQ), -1),
901  max((int)round(metaPrimitiveIt.tdc7 - shift_back * LHC_CLK_FREQ), -1),
902  max((int)round(metaPrimitiveIt.tdc8 - shift_back * LHC_CLK_FREQ), -1)};
903 
904  int pathLat[8] = {metaPrimitiveIt.lat1,
905  metaPrimitiveIt.lat2,
906  metaPrimitiveIt.lat3,
907  metaPrimitiveIt.lat4,
908  metaPrimitiveIt.lat5,
909  metaPrimitiveIt.lat6,
910  metaPrimitiveIt.lat7,
911  metaPrimitiveIt.lat8};
912 
913  // phiTP (extended DF)
914  outExtP2Ph.emplace_back(
915  L1Phase2MuDTExtPhDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
916  chId.wheel(), // uwh (m_wheel)
917  sectorTP, // usc (m_sector)
918  chId.station(), // ust (m_station)
919  sl, // ust (m_station)
920  (int)round(metaPrimitiveIt.phi * PHIRES_CONV), // uphi (m_phiAngle)
921  (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV), // uphib (m_phiBending)
922  metaPrimitiveIt.quality, // uqua (m_qualityCode)
923  metaPrimitiveIt.index, // uind (m_segmentIndex)
924  tp_t0, // ut0 (m_t0Segment)
925  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
926  (int)round(metaPrimitiveIt.x * 1000), // ux (m_xLocal)
927  (int)round(metaPrimitiveIt.tanPhi * 1000), // utan (m_tanPsi)
928  (int)round(metaPrimitiveIt.phi_cmssw * PHIRES_CONV), // uphi (m_phiAngleCMSSW)
929  (int)round(metaPrimitiveIt.phiB_cmssw * PHIBRES_CONV), // uphib (m_phiBendingCMSSW)
930  metaPrimitiveIt.rpcFlag, // urpc (m_rpcFlag)
931  pathWireId,
932  pathTDC,
933  pathLat));
934  }
935  if (df_extended_ == 0 || df_extended_ == 2) {
936  // phiTP (standard DF)
937  outP2Ph.push_back(L1Phase2MuDTPhDigi(
938  (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
939  chId.wheel(), // uwh (m_wheel)
940  sectorTP, // usc (m_sector)
941  chId.station(), // ust (m_station)
942  sl, // ust (m_station)
943  (int)round(metaPrimitiveIt.phi * PHIRES_CONV), // uphi (_phiAngle)
944  (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV), // uphib (m_phiBending)
945  metaPrimitiveIt.quality, // uqua (m_qualityCode)
946  metaPrimitiveIt.index, // uind (m_segmentIndex)
947  tp_t0, // ut0 (m_t0Segment)
948  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
949  metaPrimitiveIt.rpcFlag // urpc (m_rpcFlag)
950  ));
951  }
952  } else {
953  if (df_extended_ == 1 || df_extended_ == 2) {
954  int pathWireId[4] = {metaPrimitiveIt.wi1, metaPrimitiveIt.wi2, metaPrimitiveIt.wi3, metaPrimitiveIt.wi4};
955 
956  int pathTDC[4] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
957  max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
958  max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
959  max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1)};
960 
961  int pathLat[4] = {metaPrimitiveIt.lat1, metaPrimitiveIt.lat2, metaPrimitiveIt.lat3, metaPrimitiveIt.lat4};
962 
963  // thTP (extended DF)
964  outExtP2Th.emplace_back(
965  L1Phase2MuDTExtThDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
966  chId.wheel(), // uwh (m_wheel)
967  sectorTP, // usc (m_sector)
968  chId.station(), // ust (m_station)
969  (int)round(metaPrimitiveIt.phi * ZRES_CONV), // uz (m_zGlobal)
970  (int)round(metaPrimitiveIt.phiB * KRES_CONV), // uk (m_kSlope)
971  metaPrimitiveIt.quality, // uqua (m_qualityCode)
972  metaPrimitiveIt.index, // uind (m_segmentIndex)
973  tp_t0, // ut0 (m_t0Segment)
974  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
975  (int)round(metaPrimitiveIt.x * 1000), // ux (m_yLocal)
976  (int)round(metaPrimitiveIt.phi_cmssw * ZRES_CONV), // uphi (m_zCMSSW)
977  (int)round(metaPrimitiveIt.phiB_cmssw * KRES_CONV), // uphib (m_kCMSSW)
978  metaPrimitiveIt.rpcFlag, // urpc (m_rpcFlag)
979  pathWireId,
980  pathTDC,
981  pathLat));
982  }
983  if (df_extended_ == 0 || df_extended_ == 2) {
984  // thTP (standard DF)
985  outP2Th.push_back(L1Phase2MuDTThDigi(
986  (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
987  chId.wheel(), // uwh (m_wheel)
988  sectorTP, // usc (m_sector)
989  chId.station(), // ust (m_station)
990  (int)round(metaPrimitiveIt.phi * ZRES_CONV), // uz (m_zGlobal)
991  (int)round(metaPrimitiveIt.phiB * KRES_CONV), // uk (m_kSlope)
992  metaPrimitiveIt.quality, // uqua (m_qualityCode)
993  metaPrimitiveIt.index, // uind (m_segmentIndex)
994  tp_t0, // ut0 (m_t0Segment)
995  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
996  metaPrimitiveIt.rpcFlag // urpc (m_rpcFlag)
997  ));
998  }
999  }
1000  }
1001  }
1002 
1003  // Storing RPC hits that were not used elsewhere
1004  if (useRPC_) {
1005  for (auto rpc_dt_digi = rpc_integrator_->rpcRecHits_translated_.begin();
1006  rpc_dt_digi != rpc_integrator_->rpcRecHits_translated_.end();
1007  rpc_dt_digi++) {
1008  outP2Ph.push_back(*rpc_dt_digi);
1009  }
1010  }
1011 
1012  // Storing Phi results
1013  if (df_extended_ == 1 || df_extended_ == 2) {
1014  std::unique_ptr<L1Phase2MuDTExtPhContainer> resultExtP2Ph(new L1Phase2MuDTExtPhContainer);
1015  resultExtP2Ph->setContainer(outExtP2Ph);
1016  iEvent.put(std::move(resultExtP2Ph));
1017  }
1018  if (df_extended_ == 0 || df_extended_ == 2) {
1019  std::unique_ptr<L1Phase2MuDTPhContainer> resultP2Ph(new L1Phase2MuDTPhContainer);
1020  resultP2Ph->setContainer(outP2Ph);
1021  iEvent.put(std::move(resultP2Ph));
1022  }
1023  outExtP2Ph.clear();
1024  outExtP2Ph.erase(outExtP2Ph.begin(), outExtP2Ph.end());
1025  outP2Ph.clear();
1026  outP2Ph.erase(outP2Ph.begin(), outP2Ph.end());
1027 
1028  // Storing Theta results
1029  if (df_extended_ == 1 || df_extended_ == 2) {
1030  std::unique_ptr<L1Phase2MuDTExtThContainer> resultExtP2Th(new L1Phase2MuDTExtThContainer);
1031  resultExtP2Th->setContainer(outExtP2Th);
1032  iEvent.put(std::move(resultExtP2Th));
1033  }
1034  if (df_extended_ == 0 || df_extended_ == 2) {
1035  std::unique_ptr<L1Phase2MuDTThContainer> resultP2Th(new L1Phase2MuDTThContainer);
1036  resultP2Th->setContainer(outP2Th);
1037  iEvent.put(std::move(resultP2Th));
1038  }
1039  outExtP2Th.clear();
1040  outExtP2Th.erase(outExtP2Th.begin(), outExtP2Th.end());
1041  outP2Th.clear();
1042  outP2Th.erase(outP2Th.begin(), outP2Th.end());
1043 }
1044 
1045 void DTTrigPhase2Prod::endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
1046  grouping_obj_->finish();
1047  mpathanalyzer_->finish();
1048  mpathqualityenhancer_->finish();
1049  mpathqualityenhancerbayes_->finish();
1050  mpathredundantfilter_->finish();
1051  mpathhitsfilter_->finish();
1052  mpathassociator_->finish();
1053  rpc_integrator_->finish();
1054 };
1055 
1057  int counter = (mp.wi5 != -1) + (mp.wi6 != -1) + (mp.wi7 != -1) + (mp.wi8 != -1);
1058  return (counter > 2);
1059 }
1060 
1062  int counter = (mp.wi1 != -1) + (mp.wi2 != -1) + (mp.wi3 != -1) + (mp.wi4 != -1);
1063  return (counter > 2);
1064 }
1065 
1066 bool DTTrigPhase2Prod::hasPosRF(int wh, int sec) const { return wh > 0 || (wh == 0 && sec % 4 > 1); }
1067 
1068 void DTTrigPhase2Prod::printmP(const string& ss, const metaPrimitive& mP) const {
1069  DTSuperLayerId slId(mP.rawId);
1070  LogInfo("DTTrigPhase2Prod") << ss << (int)slId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
1071  << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
1072  << setw(5) << left << mP.tdc1 << " " << setw(5) << left << mP.tdc2 << " " << setw(5)
1073  << left << mP.tdc3 << " " << setw(5) << left << mP.tdc4 << " " << setw(10) << right
1074  << mP.x << " " << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " "
1075  << setw(13) << left << mP.chi2 << " r:" << rango(mP);
1076 }
1077 
1079  DTSuperLayerId slId(mP.rawId);
1080  LogInfo("DTTrigPhase2Prod") << (int)slId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left << mP.wi2
1081  << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " " << setw(5)
1082  << left << mP.tdc1 << " " << setw(5) << left << mP.tdc2 << " " << setw(5) << left
1083  << mP.tdc3 << " " << setw(5) << left << mP.tdc4 << " " << setw(10) << right << mP.x << " "
1084  << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13)
1085  << left << mP.chi2 << " r:" << rango(mP) << std::endl;
1086 }
1087 
1088 void DTTrigPhase2Prod::printmPC(const string& ss, const metaPrimitive& mP) const {
1089  DTChamberId ChId(mP.rawId);
1090  LogInfo("DTTrigPhase2Prod") << ss << (int)ChId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
1091  << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
1092  << setw(2) << left << mP.wi5 << " " << setw(2) << left << mP.wi6 << " " << setw(2) << left
1093  << mP.wi7 << " " << setw(2) << left << mP.wi8 << " " << setw(5) << left << mP.tdc1 << " "
1094  << setw(5) << left << mP.tdc2 << " " << setw(5) << left << mP.tdc3 << " " << setw(5)
1095  << left << mP.tdc4 << " " << setw(5) << left << mP.tdc5 << " " << setw(5) << left
1096  << mP.tdc6 << " " << setw(5) << left << mP.tdc7 << " " << setw(5) << left << mP.tdc8
1097  << " " << setw(2) << left << mP.lat1 << " " << setw(2) << left << mP.lat2 << " "
1098  << setw(2) << left << mP.lat3 << " " << setw(2) << left << mP.lat4 << " " << setw(2)
1099  << left << mP.lat5 << " " << setw(2) << left << mP.lat6 << " " << setw(2) << left
1100  << mP.lat7 << " " << setw(2) << left << mP.lat8 << " " << setw(10) << right << mP.x << " "
1101  << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13)
1102  << left << mP.chi2 << " r:" << rango(mP);
1103 }
1104 
1106  DTChamberId ChId(mP.rawId);
1107  LogInfo("DTTrigPhase2Prod") << (int)ChId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left << mP.wi2
1108  << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " " << setw(2)
1109  << left << mP.wi5 << " " << setw(2) << left << mP.wi6 << " " << setw(2) << left << mP.wi7
1110  << " " << setw(2) << left << mP.wi8 << " " << setw(5) << left << mP.tdc1 << " " << setw(5)
1111  << left << mP.tdc2 << " " << setw(5) << left << mP.tdc3 << " " << setw(5) << left
1112  << mP.tdc4 << " " << setw(5) << left << mP.tdc5 << " " << setw(5) << left << mP.tdc6
1113  << " " << setw(5) << left << mP.tdc7 << " " << setw(5) << left << mP.tdc8 << " "
1114  << setw(2) << left << mP.lat1 << " " << setw(2) << left << mP.lat2 << " " << setw(2)
1115  << left << mP.lat3 << " " << setw(2) << left << mP.lat4 << " " << setw(2) << left
1116  << mP.lat5 << " " << setw(2) << left << mP.lat6 << " " << setw(2) << left << mP.lat7
1117  << " " << setw(2) << left << mP.lat8 << " " << setw(10) << right << mP.x << " " << setw(9)
1118  << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13) << left
1119  << mP.chi2 << " r:" << rango(mP) << std::endl;
1120 }
1121 
1123  if (mp.quality == 1 or mp.quality == 2)
1124  return 3;
1125  if (mp.quality == 3 or mp.quality == 4)
1126  return 4;
1127  return mp.quality;
1128 }
1129 
1130 void DTTrigPhase2Prod::assignIndex(std::vector<metaPrimitive>& inMPaths) {
1131  std::map<int, std::vector<metaPrimitive>> primsPerBX;
1132  for (const auto& metaPrimitive : inMPaths) {
1133  int BX = round(metaPrimitive.t0 / 25.);
1134  primsPerBX[BX].push_back(metaPrimitive);
1135  }
1136  inMPaths.clear();
1137  for (auto& prims : primsPerBX) {
1138  assignIndexPerBX(prims.second);
1139  for (const auto& primitive : prims.second)
1140  if (primitive.index <= max_index_)
1141  inMPaths.push_back(primitive);
1142  }
1143 }
1144 
1145 void DTTrigPhase2Prod::assignIndexPerBX(std::vector<metaPrimitive>& inMPaths) {
1146  // First we asociate a new index to the metaprimitive depending on quality or phiB;
1147  uint32_t rawId = -1;
1148  int numP = -1;
1149  for (auto& metaPrimitiveIt : inMPaths) {
1150  numP++;
1151  rawId = metaPrimitiveIt.rawId;
1152  int iOrder = assignQualityOrder(metaPrimitiveIt);
1153  int inf = 0;
1154  int numP2 = -1;
1155  for (auto& metaPrimitiveItN : inMPaths) {
1156  int nOrder = assignQualityOrder(metaPrimitiveItN);
1157  numP2++;
1158  if (rawId != metaPrimitiveItN.rawId)
1159  continue;
1160  if (numP2 == numP) {
1161  metaPrimitiveIt.index = inf;
1162  break;
1163  } else if (iOrder < nOrder) {
1164  inf++;
1165  } else if (iOrder > nOrder) {
1166  metaPrimitiveItN.index++;
1167  } else if (iOrder == nOrder) {
1168  if (std::abs(metaPrimitiveIt.phiB) >= std::abs(metaPrimitiveItN.phiB)) {
1169  inf++;
1170  } else if (std::abs(metaPrimitiveIt.phiB) < std::abs(metaPrimitiveItN.phiB)) {
1171  metaPrimitiveItN.index++;
1172  }
1173  }
1174  } // ending second for
1175  } // ending first for
1176 }
1177 
1179  if (mP.quality > 8 || mP.quality < 1)
1180  return -1;
1181 
1182  return qmap_.find(mP.quality)->second;
1183 }
1184 
1185 std::vector<DTDigiCollection*> DTTrigPhase2Prod::distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ) {
1186  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*> tmpVector;
1187  tmpVector.clear();
1188  std::vector<DTDigiCollection*> collVector;
1189  collVector.clear();
1190  while (!inQ.empty()) {
1191  processDigi(inQ, tmpVector);
1192  }
1193  for (auto& sQ : tmpVector) {
1194  DTDigiCollection tmpColl;
1195  while (!sQ->empty()) {
1196  tmpColl.insertDigi((sQ->front().first), (sQ->front().second));
1197  sQ->pop();
1198  }
1199  collVector.push_back(&tmpColl);
1200  }
1201  return collVector;
1202 }
1203 
1204 void DTTrigPhase2Prod::processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
1205  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec) {
1206  bool classified = false;
1207  if (!vec.empty()) {
1208  for (auto& sC : vec) { // Conditions for entering a super cell.
1209  if ((sC->front().second.time() + superCelltimewidth_) > inQ.front().second.time()) {
1210  // Time requirement
1211  if (TMath::Abs(sC->front().second.wire() - inQ.front().second.wire()) <= superCellhalfspacewidth_) {
1212  // Spatial requirement
1213  sC->push(std::move(inQ.front()));
1214  classified = true;
1215  }
1216  }
1217  }
1218  }
1219  if (classified) {
1220  inQ.pop();
1221  return;
1222  }
1223 
1224  std::queue<std::pair<DTLayerId, DTDigi>> newQueue;
1225 
1226  std::pair<DTLayerId, DTDigi> tmpPair;
1227  tmpPair = std::move(inQ.front());
1228  newQueue.push(tmpPair);
1229  inQ.pop();
1230 
1231  vec.push_back(&newQueue);
1232 }
1233 
1235  // dtTriggerPhase2PrimitiveDigis
1237  desc.add<edm::InputTag>("digiTag", edm::InputTag("CalibratedDigis"));
1238  desc.add<int>("timeTolerance", 999999);
1239  desc.add<double>("tanPhiTh", 1.0);
1240  desc.add<double>("tanPhiThw2max", 1.3);
1241  desc.add<double>("tanPhiThw2min", 0.5);
1242  desc.add<double>("tanPhiThw1max", 0.9);
1243  desc.add<double>("tanPhiThw1min", 0.2);
1244  desc.add<double>("tanPhiThw0", 0.5);
1245  desc.add<double>("chi2Th", 0.01);
1246  desc.add<double>("chi2corTh", 0.1);
1247  desc.add<bool>("useBX_correlation", false);
1248  desc.add<double>("dT0_correlate_TP", 25.0);
1249  desc.add<int>("dBX_correlate_TP", 0);
1250  desc.add<double>("dTanPsi_correlate_TP", 99999.0);
1251  desc.add<bool>("clean_chi2_correlation", true);
1252  desc.add<bool>("allow_confirmation", true);
1253  desc.add<double>("minx_match_2digis", 1.0);
1254  desc.add<int>("scenario", 0);
1255  desc.add<int>("df_extended", 0);
1256  desc.add<int>("max_primitives", 999);
1257  desc.add<bool>("output_mixer", false);
1258  desc.add<bool>("output_latpredictor", false);
1259  desc.add<bool>("output_slfitter", false);
1260  desc.add<bool>("output_slfilter", false);
1261  desc.add<bool>("output_confirmed", false);
1262  desc.add<bool>("output_matcher", false);
1263  desc.add<edm::FileInPath>("ttrig_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_ttrig.txt"));
1264  desc.add<edm::FileInPath>("z_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_z.txt"));
1265  desc.add<edm::FileInPath>("lut_sl1", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_sl1.dat"));
1266  desc.add<edm::FileInPath>("lut_sl2", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_slx.dat"));
1267  desc.add<edm::FileInPath>("lut_sl3", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_sl3.dat"));
1268  desc.add<edm::FileInPath>("lut_2sl", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_2sl.dat"));
1269  desc.add<edm::FileInPath>("shift_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_x.txt"));
1270  desc.add<edm::FileInPath>("maxdrift_filename",
1271  edm::FileInPath("L1Trigger/DTTriggerPhase2/data/drift_time_per_chamber.txt"));
1272  desc.add<edm::FileInPath>("shift_theta_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/theta_shift.txt"));
1273  desc.add<edm::FileInPath>("global_coords_filename",
1274  edm::FileInPath("L1Trigger/DTTriggerPhase2/data/global_coord_perp_x_phi0.txt"));
1275  desc.add<edm::FileInPath>("laterality_filename",
1276  edm::FileInPath("L1Trigger/DTTriggerPhase2/data/lat_predictions.dat"));
1277  desc.add<int>("algo", 0);
1278  desc.add<int>("minHits4Fit", 3);
1279  desc.add<bool>("splitPathPerSL", true);
1280  desc.addUntracked<bool>("debug", false);
1281  desc.addUntracked<bool>("dump", false);
1282  desc.add<edm::InputTag>("rpcRecHits", edm::InputTag("rpcRecHits"));
1283  desc.add<bool>("useRPC", false);
1284  desc.add<int>("bx_window", 1);
1285  desc.add<double>("phi_window", 50.0);
1286  desc.add<int>("max_quality_to_overwrite_t0", 9);
1287  desc.add<bool>("storeAllRPCHits", false);
1288  desc.add<bool>("activateBuffer", false);
1289  desc.add<double>("superCelltimewidth", 400);
1290  desc.add<int>("superCellspacewidth", 20);
1291  {
1293  psd0.addUntracked<bool>("debug", false);
1294  psd0.add<double>("angletan", 0.3);
1295  psd0.add<double>("anglebinwidth", 1.0);
1296  psd0.add<double>("posbinwidth", 2.1);
1297  psd0.add<double>("maxdeltaAngDeg", 10);
1298  psd0.add<double>("maxdeltaPos", 10);
1299  psd0.add<int>("UpperNumber", 6);
1300  psd0.add<int>("LowerNumber", 4);
1301  psd0.add<double>("MaxDistanceToWire", 0.03);
1302  psd0.add<int>("minNLayerHits", 6);
1303  psd0.add<int>("minSingleSLHitsMax", 3);
1304  psd0.add<int>("minSingleSLHitsMin", 3);
1305  psd0.add<bool>("allowUncorrelatedPatterns", true);
1306  psd0.add<int>("minUncorrelatedHits", 3);
1307  desc.add<edm::ParameterSetDescription>("HoughGrouping", psd0);
1308  }
1309  {
1311  psd0.add<edm::FileInPath>(
1312  "pattern_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/PseudoBayesPatterns_uncorrelated_v0.root"));
1313  psd0.addUntracked<bool>("debug", false);
1314  psd0.add<int>("minNLayerHits", 3);
1315  psd0.add<int>("minSingleSLHitsMax", 3);
1316  psd0.add<int>("minSingleSLHitsMin", 0);
1317  psd0.add<int>("allowedVariance", 1);
1318  psd0.add<bool>("allowDuplicates", false);
1319  psd0.add<bool>("setLateralities", true);
1320  psd0.add<bool>("allowUncorrelatedPatterns", true);
1321  psd0.add<int>("minUncorrelatedHits", 3);
1322  psd0.add<bool>("saveOnPlace", true);
1323  psd0.add<int>("maxPathsPerMatch", 256);
1324  desc.add<edm::ParameterSetDescription>("PseudoBayesPattern", psd0);
1325  }
1326  descriptions.add("dtTriggerPhase2PrimitiveDigis", desc);
1327 }
1328 
int station() const
Return the station number.
Definition: DTChamberId.h:45
std::unique_ptr< MotherGrouping > grouping_obj_
std::map< DTChamberId, DTDigiCollection, std::less< DTChamberId > > DTDigiMap
bool hasPosRF(int wh, int sec) const
void produce(edm::Event &iEvent, const edm::EventSetup &iEventSetup) override
Producer: process every event and generates trigger data.
int superLayer() const
Return the superlayer number.
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
constexpr int CHI2RES_CONV
Definition: constants.h:365
constexpr float PHIBRES_CONV
Definition: constants.h:364
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:32
std::unique_ptr< LateralityProvider > latprovider_
DTTrigPhase2Prod(const edm::ParameterSet &pset)
Constructor.
void processDigi(std::queue< std::pair< DTLayerId, DTDigi >> &inQ, std::vector< std::queue< std::pair< DTLayerId, DTDigi >> *> &vec)
std::unique_ptr< MPFilter > mpathhitsfilter_
edm::EDGetTokenT< DTDigiCollection > dtDigisToken_
std::vector< DTDigiCollection * > distribDigis(std::queue< std::pair< DTLayerId, DTDigi >> &inQ)
std::unique_ptr< RPCIntegrator > rpc_integrator_
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomH
std::unique_ptr< MuonPathConfirmator > mpathconfirmator_
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsLabel_
Primitive< F, X >::type primitive(const F &f)
Definition: Primitive.h:41
std::vector< short > latcomb
void assignIndex(std::vector< metaPrimitive > &inMPaths)
const DTGeometry * dtGeo_
int iEvent
Definition: GenABIO.cc:224
void printmPC(const std::string &ss, const metaPrimitive &mP) const
DTChamberId chamberId() const
Return the corresponding ChamberId.
bool inner(const metaPrimitive &mp) const
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::unique_ptr< MuonPathAnalyzer > mpathanalyzer_
constexpr float KRES_CONV
Definition: constants.h:368
std::unique_ptr< MPFilter > mpathqualityenhancer_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void printmP(const std::string &ss, const metaPrimitive &mP) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr float ZRES_CONV
Definition: constants.h:367
ParameterDescriptionBase * add(U const &iLabel, T const &value)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::unique_ptr< MPFilter > mpathcorfilter_
RunID const & id() const
Definition: RunBase.h:39
Log< level::Info, false > LogInfo
MP_QUALITY
Definition: constants.h:44
bool outer(const metaPrimitive &mp) const
void endRun(edm::Run const &iRun, const edm::EventSetup &iEventSetup) override
endRun: finish things
std::unique_ptr< MPFilter > mpathredundantfilter_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr int TIME_TO_TDC_COUNTS
Definition: constants.h:235
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::shared_ptr< GlobalCoordsObtainer > globalcoordsobtainer_
std::pair< const_iterator, const_iterator > Range
std::unique_ptr< MuonPathAnalyzer > mpathassociator_
double b
Definition: hdecay.h:120
std::vector< DigiType >::const_iterator const_iterator
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::unique_ptr< MPFilter > mpathqualityenhancerbayes_
constexpr int LHC_CLK_FREQ
Definition: constants.h:222
int rango(const metaPrimitive &mp) const
DTDigiMap::iterator DTDigiMap_iterator
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
HLT enums.
int sector() const
Definition: DTChamberId.h:52
void assignIndexPerBX(std::vector< metaPrimitive > &inMPaths)
~DTTrigPhase2Prod() override
Destructor.
double a
Definition: hdecay.h:121
DTDigiMap::const_iterator DTDigiMap_const_iterator
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:48
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
constexpr float PHIRES_CONV
Definition: constants.h:363
int assignQualityOrder(const metaPrimitive &mP) const
const std::unordered_map< int, int > qmap_
void beginRun(edm::Run const &iRun, const edm::EventSetup &iEventSetup) override
Create Trigger Units before starting event processing.
def move(src, dest)
Definition: eostools.py:511
std::vector< std::pair< int, MuonPath > > primitives_
RunNumber_t run() const
Definition: RunID.h:26
Definition: Run.h:45
#define LogDebug(id)