CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DTTrigPhase2Prod.cc
Go to the documentation of this file.
10 
17 
21 
24 
37 
51 
52 // DT trigger GeomUtils
54 
55 //RPC TP
60 
61 #include <fstream>
62 #include <iostream>
63 #include <queue>
64 #include <cmath>
65 
66 using namespace edm;
67 using namespace std;
68 using namespace cmsdt;
69 
71  typedef std::map<DTChamberId, DTDigiCollection, std::less<DTChamberId>> DTDigiMap;
72  typedef DTDigiMap::iterator DTDigiMap_iterator;
73  typedef DTDigiMap::const_iterator DTDigiMap_const_iterator;
74 
75 public:
78 
80  ~DTTrigPhase2Prod() override;
81 
83  void beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
84 
86  void produce(edm::Event& iEvent, const edm::EventSetup& iEventSetup) override;
87 
89  void endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
90 
91  // Methods
92  int rango(const metaPrimitive& mp) const;
93  bool outer(const metaPrimitive& mp) const;
94  bool inner(const metaPrimitive& mp) const;
95  void printmP(const std::string& ss, const metaPrimitive& mP) const;
96  void printmPC(const std::string& ss, const metaPrimitive& mP) const;
97  bool hasPosRF(int wh, int sec) const;
98 
99  // Getter-methods
100  MP_QUALITY getMinimumQuality(void);
101 
102  // Setter-methods
103  void setChiSquareThreshold(float ch2Thr);
104  void setMinimumQuality(MP_QUALITY q);
105 
106  // data-members
109  std::vector<std::pair<int, MuonPath>> primitives_;
110 
111 private:
112  // Trigger Configuration Manager CCB validity flag
114 
115  // BX offset used to correct DTTPG output
117 
118  // Debug Flag
119  bool debug_;
120  bool dump_;
126 
127  // ParameterSet
130 
131  // Grouping attributes and methods
132  int algo_; // Grouping code
133  std::unique_ptr<MotherGrouping> grouping_obj_;
134  std::unique_ptr<MuonPathAnalyzer> mpathanalyzer_;
135  std::unique_ptr<MPFilter> mpathqualityenhancer_;
136  std::unique_ptr<MPFilter> mpathredundantfilter_;
137  // std::unique_ptr<MPFilter> mpathhitsfilter_;
138  std::unique_ptr<MuonPathAssociator> mpathassociator_;
139  std::shared_ptr<GlobalCoordsObtainer> globalcoordsobtainer_;
140 
141  // Buffering
145  std::vector<DTDigiCollection*> distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ);
146  void processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
147  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec);
148 
149  // RPC
150  std::unique_ptr<RPCIntegrator> rpc_integrator_;
151  bool useRPC_;
152 
153  void assignIndex(std::vector<metaPrimitive>& inMPaths);
154  void assignIndexPerBX(std::vector<metaPrimitive>& inMPaths);
155  int assignQualityOrder(const metaPrimitive& mP) const;
156 
157  const std::unordered_map<int, int> qmap_;
158 };
159 
160 namespace {
161  struct {
162  bool operator()(std::pair<DTLayerId, DTDigi> a, std::pair<DTLayerId, DTDigi> b) const {
163  return (a.second.time() < b.second.time());
164  }
165  } const DigiTimeOrdering;
166 } // namespace
167 
169  : qmap_({{9, 9}, {8, 8}, {7, 6}, {6, 7}, {5, 3}, {4, 5}, {3, 4}, {2, 2}, {1, 1}}) {
170  produces<L1Phase2MuDTPhContainer>();
174 
175  debug_ = pset.getUntrackedParameter<bool>("debug");
176  dump_ = pset.getUntrackedParameter<bool>("dump");
177 
178  do_correlation_ = pset.getParameter<bool>("do_correlation");
179  scenario_ = pset.getParameter<int>("scenario");
180 
181  df_extended_ = pset.getParameter<int>("df_extended");
182 
183  dtDigisToken_ = consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiTag"));
184 
185  rpcRecHitsLabel_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("rpcRecHits"));
186  useRPC_ = pset.getParameter<bool>("useRPC");
187 
188  // Choosing grouping scheme:
189  algo_ = pset.getParameter<int>("algo");
190 
191  // Local to global coordinates approach
192  geometry_tag_ = pset.getUntrackedParameter<std::string>("geometry_tag", "");
193 
194  edm::ConsumesCollector consumesColl(consumesCollector());
195  globalcoordsobtainer_ = std::make_shared<GlobalCoordsObtainer>(pset);
196  globalcoordsobtainer_->generate_luts();
197 
198  if (algo_ == PseudoBayes) {
199  grouping_obj_ =
200  std::make_unique<PseudoBayesGrouping>(pset.getParameter<edm::ParameterSet>("PseudoBayesPattern"), consumesColl);
201  } else if (algo_ == HoughTrans) {
202  grouping_obj_ =
203  std::make_unique<HoughGrouping>(pset.getParameter<edm::ParameterSet>("HoughGrouping"), consumesColl);
204  } else {
205  grouping_obj_ = std::make_unique<InitialGrouping>(pset, consumesColl);
206  }
207 
208  if (algo_ == Standard) {
209  if (debug_)
210  LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: JM analyzer";
211  mpathanalyzer_ = std::make_unique<MuonPathAnalyticAnalyzer>(pset, consumesColl, globalcoordsobtainer_);
212  } else {
213  if (debug_)
214  LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: Full chamber analyzer";
215  mpathanalyzer_ = std::make_unique<MuonPathAnalyzerInChamber>(pset, consumesColl, globalcoordsobtainer_);
216  }
217 
218  // Getting buffer option
219  activateBuffer_ = pset.getParameter<bool>("activateBuffer");
220  superCellhalfspacewidth_ = pset.getParameter<int>("superCellspacewidth") / 2;
221  superCelltimewidth_ = pset.getParameter<double>("superCelltimewidth");
222 
223  mpathqualityenhancer_ = std::make_unique<MPQualityEnhancerFilter>(pset);
224  mpathredundantfilter_ = std::make_unique<MPRedundantFilter>(pset);
225  mpathassociator_ = std::make_unique<MuonPathAssociator>(pset, consumesColl, globalcoordsobtainer_);
226  rpc_integrator_ = std::make_unique<RPCIntegrator>(pset, consumesColl);
227 
228  dtGeomH = esConsumes<DTGeometry, MuonGeometryRecord, edm::Transition::BeginRun>();
229 }
230 
232  if (debug_)
233  LogDebug("DTTrigPhase2Prod") << "DTp2: calling destructor" << std::endl;
234 }
235 
236 void DTTrigPhase2Prod::beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
237  if (debug_)
238  LogDebug("DTTrigPhase2Prod") << "beginRun " << iRun.id().run();
239  if (debug_)
240  LogDebug("DTTrigPhase2Prod") << "beginRun: getting DT geometry";
241 
242  grouping_obj_->initialise(iEventSetup); // Grouping object initialisation
243  mpathanalyzer_->initialise(iEventSetup); // Analyzer object initialisation
244  mpathqualityenhancer_->initialise(iEventSetup); // Filter object initialisation
245  mpathredundantfilter_->initialise(iEventSetup); // Filter object initialisation
246  mpathassociator_->initialise(iEventSetup); // Associator object initialisation
247 
249  iEventSetup.get<MuonGeometryRecord>().get(geometry_tag_, geom);
250  dtGeo_ = &(*geom);
251 }
252 
253 void DTTrigPhase2Prod::produce(Event& iEvent, const EventSetup& iEventSetup) {
254  if (debug_)
255  LogDebug("DTTrigPhase2Prod") << "produce";
257  iEvent.getByToken(dtDigisToken_, dtdigis);
258 
259  if (debug_)
260  LogDebug("DTTrigPhase2Prod") << "\t Getting the RPC RecHits" << std::endl;
262  iEvent.getByToken(rpcRecHitsLabel_, rpcRecHits);
263 
265  // GROUPING CODE:
267 
268  DTDigiMap digiMap;
270  for (const auto& detUnitIt : *dtdigis) {
271  const DTLayerId& layId = detUnitIt.first;
272  const DTChamberId chambId = layId.superlayerId().chamberId();
273  const DTDigiCollection::Range& range = detUnitIt.second;
274  digiMap[chambId].put(range, layId);
275  }
276 
277  // generate a list muon paths for each event!!!
278  if (debug_ && activateBuffer_)
279  LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber using a buffer and super cells.";
280  else if (debug_)
281  LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber.";
282 
283  MuonPathPtrs muonpaths;
284  for (const auto& ich : dtGeo_->chambers()) {
285  // 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.
286  const DTChamber* chamb = ich;
287  DTChamberId chid = chamb->id();
288  DTDigiMap_iterator dmit = digiMap.find(chid);
289 
290  if (dmit == digiMap.end())
291  continue;
292 
293  if (activateBuffer_) { // Use buffering (per chamber) or not
294  // Import digis from the station
295  std::vector<std::pair<DTLayerId, DTDigi>> tmpvec;
296  tmpvec.clear();
297 
298  for (const auto& dtLayerIdIt : (*dmit).second) {
299  for (DTDigiCollection::const_iterator digiIt = (dtLayerIdIt.second).first;
300  digiIt != (dtLayerIdIt.second).second;
301  digiIt++) {
302  tmpvec.emplace_back(dtLayerIdIt.first, *digiIt);
303  }
304  }
305 
306  // Check to enhance CPU time usage
307  if (tmpvec.empty())
308  continue;
309 
310  // Order digis depending on TDC time and insert them into a queue (FIFO buffer). TODO: adapt for MC simulations.
311  std::sort(tmpvec.begin(), tmpvec.end(), DigiTimeOrdering);
312  std::queue<std::pair<DTLayerId, DTDigi>> timequeue;
313 
314  for (const auto& elem : tmpvec)
315  timequeue.emplace(elem);
316  tmpvec.clear();
317 
318  // Distribute the digis from the queue into supercells
319  std::vector<DTDigiCollection*> superCells;
320  superCells = distribDigis(timequeue);
321 
322  // Process each supercell & collect the resulting muonpaths (as the muonpaths std::vector is only enlarged each time
323  // the groupings access it, it's not needed to "collect" the final products).
324 
325  while (!superCells.empty()) {
326  grouping_obj_->run(iEvent, iEventSetup, *(superCells.back()), muonpaths);
327  superCells.pop_back();
328  }
329  } else {
330  grouping_obj_->run(iEvent, iEventSetup, (*dmit).second, muonpaths);
331  }
332  }
333  digiMap.clear();
334 
335  if (dump_) {
336  for (unsigned int i = 0; i < muonpaths.size(); i++) {
337  stringstream ss;
338  ss << iEvent.id().event() << " mpath " << i << ": ";
339  for (int lay = 0; lay < muonpaths.at(i)->nprimitives(); lay++)
340  ss << muonpaths.at(i)->primitive(lay)->channelId() << " ";
341  for (int lay = 0; lay < muonpaths.at(i)->nprimitives(); lay++)
342  ss << muonpaths.at(i)->primitive(lay)->tdcTimeStamp() << " ";
343  for (int lay = 0; lay < muonpaths.at(i)->nprimitives(); lay++)
344  ss << muonpaths.at(i)->primitive(lay)->laterality() << " ";
345  LogInfo("DTTrigPhase2Prod") << ss.str();
346  }
347  }
348 
349  // FILTER GROUPING
350  MuonPathPtrs filteredmuonpaths;
351  if (algo_ == Standard) {
352  mpathredundantfilter_->run(iEvent, iEventSetup, muonpaths, filteredmuonpaths);
353  }
354 
355  if (dump_) {
356  for (unsigned int i = 0; i < filteredmuonpaths.size(); i++) {
357  stringstream ss;
358  ss << iEvent.id().event() << " filt. mpath " << i << ": ";
359  for (int lay = 0; lay < filteredmuonpaths.at(i)->nprimitives(); lay++)
360  ss << filteredmuonpaths.at(i)->primitive(lay)->channelId() << " ";
361  for (int lay = 0; lay < filteredmuonpaths.at(i)->nprimitives(); lay++)
362  ss << filteredmuonpaths.at(i)->primitive(lay)->tdcTimeStamp() << " ";
363  LogInfo("DTTrigPhase2Prod") << ss.str();
364  }
365  }
366 
370 
371  if (debug_)
372  LogDebug("DTTrigPhase2Prod") << "MUON PATHS found: " << muonpaths.size() << " (" << filteredmuonpaths.size()
373  << ") in event " << iEvent.id().event();
374  if (debug_)
375  LogDebug("DTTrigPhase2Prod") << "filling NmetaPrimtives" << std::endl;
376  std::vector<metaPrimitive> metaPrimitives;
377  MuonPathPtrs outmpaths;
378  if (algo_ == Standard) {
379  if (debug_)
380  LogDebug("DTTrigPhase2Prod") << "Fitting 1SL ";
381  mpathanalyzer_->run(iEvent, iEventSetup, filteredmuonpaths, metaPrimitives);
382  } else {
383  // implementation for advanced (2SL) grouping, no filter required..
384  if (debug_)
385  LogDebug("DTTrigPhase2Prod") << "Fitting 2SL at once ";
386  mpathanalyzer_->run(iEvent, iEventSetup, muonpaths, outmpaths);
387  }
388 
389  if (dump_) {
390  for (unsigned int i = 0; i < outmpaths.size(); i++) {
391  LogInfo("DTTrigPhase2Prod") << iEvent.id().event() << " mp " << i << ": " << outmpaths.at(i)->bxTimeValue() << " "
392  << outmpaths.at(i)->horizPos() << " " << outmpaths.at(i)->tanPhi() << " "
393  << outmpaths.at(i)->phi() << " " << outmpaths.at(i)->phiB() << " "
394  << outmpaths.at(i)->quality() << " " << outmpaths.at(i)->chiSquare();
395  }
396  for (unsigned int i = 0; i < metaPrimitives.size(); i++) {
397  stringstream ss;
398  ss << iEvent.id().event() << " mp " << i << ": ";
399  printmP(ss.str(), metaPrimitives.at(i));
400  }
401  }
402 
403  muonpaths.clear();
404  filteredmuonpaths.clear();
405 
407  // FILTER SECTIONS:
409 
410  if (debug_)
411  LogDebug("DTTrigPhase2Prod") << "declaring new vector for filtered" << std::endl;
412 
413  std::vector<metaPrimitive> filteredMetaPrimitives;
414  if (algo_ == Standard)
415  mpathqualityenhancer_->run(iEvent, iEventSetup, metaPrimitives, filteredMetaPrimitives);
416 
417  if (dump_) {
418  for (unsigned int i = 0; i < filteredMetaPrimitives.size(); i++) {
419  stringstream ss;
420  ss << iEvent.id().event() << " filtered mp " << i << ": ";
421  printmP(ss.str(), filteredMetaPrimitives.at(i));
422  }
423  }
424 
425  metaPrimitives.clear();
426  metaPrimitives.erase(metaPrimitives.begin(), metaPrimitives.end());
427 
428  if (debug_)
429  LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
430  << filteredMetaPrimitives.size() << " filteredMetaPrimitives (superlayer)"
431  << std::endl;
432  if (debug_)
433  LogDebug("DTTrigPhase2Prod") << "filteredMetaPrimitives: starting correlations" << std::endl;
434 
438 
439  std::vector<metaPrimitive> correlatedMetaPrimitives;
440  if (algo_ == Standard)
441  mpathassociator_->run(iEvent, iEventSetup, dtdigis, filteredMetaPrimitives, correlatedMetaPrimitives);
442  else {
443  for (const auto& muonpath : outmpaths) {
444  correlatedMetaPrimitives.emplace_back(muonpath->rawId(),
445  (double)muonpath->bxTimeValue(),
446  muonpath->horizPos(),
447  muonpath->tanPhi(),
448  muonpath->phi(),
449  muonpath->phiB(),
450  muonpath->phi_cmssw(),
451  muonpath->phiB_cmssw(),
452  muonpath->chiSquare(),
453  (int)muonpath->quality(),
454  muonpath->primitive(0)->channelId(),
455  muonpath->primitive(0)->tdcTimeStamp(),
456  muonpath->primitive(0)->laterality(),
457  muonpath->primitive(1)->channelId(),
458  muonpath->primitive(1)->tdcTimeStamp(),
459  muonpath->primitive(1)->laterality(),
460  muonpath->primitive(2)->channelId(),
461  muonpath->primitive(2)->tdcTimeStamp(),
462  muonpath->primitive(2)->laterality(),
463  muonpath->primitive(3)->channelId(),
464  muonpath->primitive(3)->tdcTimeStamp(),
465  muonpath->primitive(3)->laterality(),
466  muonpath->primitive(4)->channelId(),
467  muonpath->primitive(4)->tdcTimeStamp(),
468  muonpath->primitive(4)->laterality(),
469  muonpath->primitive(5)->channelId(),
470  muonpath->primitive(5)->tdcTimeStamp(),
471  muonpath->primitive(5)->laterality(),
472  muonpath->primitive(6)->channelId(),
473  muonpath->primitive(6)->tdcTimeStamp(),
474  muonpath->primitive(6)->laterality(),
475  muonpath->primitive(7)->channelId(),
476  muonpath->primitive(7)->tdcTimeStamp(),
477  muonpath->primitive(7)->laterality());
478  }
479  }
480  filteredMetaPrimitives.clear();
481 
482  if (debug_)
483  LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
484  << correlatedMetaPrimitives.size() << " correlatedMetPrimitives (chamber)";
485 
486  if (dump_) {
487  LogInfo("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
488  << correlatedMetaPrimitives.size() << " correlatedMetPrimitives (chamber)";
489 
490  for (unsigned int i = 0; i < correlatedMetaPrimitives.size(); i++) {
491  stringstream ss;
492  ss << iEvent.id().event() << " correlated mp " << i << ": ";
493  printmPC(ss.str(), correlatedMetaPrimitives.at(i));
494  }
495  }
496 
497  double shift_back = 0;
498  if (scenario_ == MC) //scope for MC
499  shift_back = 400;
500  else if (scenario_ == DATA) //scope for data
501  shift_back = 0;
502  else if (scenario_ == SLICE_TEST) //scope for slice test
503  shift_back = 400;
504 
505  // RPC integration
506  if (useRPC_) {
507  rpc_integrator_->initialise(iEventSetup, shift_back);
508  rpc_integrator_->prepareMetaPrimitives(rpcRecHits);
509  rpc_integrator_->matchWithDTAndUseRPCTime(correlatedMetaPrimitives);
510  rpc_integrator_->makeRPCOnlySegments();
511  rpc_integrator_->storeRPCSingleHits();
512  rpc_integrator_->removeRPCHitsUsed();
513  }
514 
516  vector<L1Phase2MuDTPhDigi> outP2Ph;
517  vector<L1Phase2MuDTExtPhDigi> outExtP2Ph;
518  vector<L1Phase2MuDTThDigi> outP2Th;
519  vector<L1Phase2MuDTExtThDigi> outExtP2Th;
520 
521  // Assigning index value
522  assignIndex(correlatedMetaPrimitives);
523  for (const auto& metaPrimitiveIt : correlatedMetaPrimitives) {
524  DTChamberId chId(metaPrimitiveIt.rawId);
525  DTSuperLayerId slId(metaPrimitiveIt.rawId);
526  if (debug_)
527  LogDebug("DTTrigPhase2Prod") << "looping in final vector: SuperLayerId" << chId << " x=" << metaPrimitiveIt.x
528  << " quality=" << metaPrimitiveIt.quality
529  << " BX=" << round(metaPrimitiveIt.t0 / 25.) << " index=" << metaPrimitiveIt.index;
530 
531  int sectorTP = chId.sector();
532  //sectors 13 and 14 exist only for the outermost stations for sectors 4 and 10 respectively
533  //due to the larger MB4 that are divided into two.
534  if (sectorTP == 13)
535  sectorTP = 4;
536  if (sectorTP == 14)
537  sectorTP = 10;
538  sectorTP = sectorTP - 1;
539  int sl = 0;
540  if (metaPrimitiveIt.quality < LOWLOWQ || metaPrimitiveIt.quality == CHIGHQ) {
541  if (inner(metaPrimitiveIt))
542  sl = 1;
543  else
544  sl = 3;
545  }
546 
547  if (debug_)
548  LogDebug("DTTrigPhase2Prod") << "pushing back phase-2 dataformat carlo-federica dataformat";
549 
550  if (slId.superLayer() != 2) {
551  if (df_extended_ == 1 || df_extended_ == 2) {
552  int pathWireId[8] = {metaPrimitiveIt.wi1,
553  metaPrimitiveIt.wi2,
554  metaPrimitiveIt.wi3,
555  metaPrimitiveIt.wi4,
556  metaPrimitiveIt.wi5,
557  metaPrimitiveIt.wi6,
558  metaPrimitiveIt.wi7,
559  metaPrimitiveIt.wi8};
560 
561  int pathTDC[8] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
562  max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
563  max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
564  max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1),
565  max((int)round(metaPrimitiveIt.tdc5 - shift_back * LHC_CLK_FREQ), -1),
566  max((int)round(metaPrimitiveIt.tdc6 - shift_back * LHC_CLK_FREQ), -1),
567  max((int)round(metaPrimitiveIt.tdc7 - shift_back * LHC_CLK_FREQ), -1),
568  max((int)round(metaPrimitiveIt.tdc8 - shift_back * LHC_CLK_FREQ), -1)};
569 
570  int pathLat[8] = {metaPrimitiveIt.lat1,
571  metaPrimitiveIt.lat2,
572  metaPrimitiveIt.lat3,
573  metaPrimitiveIt.lat4,
574  metaPrimitiveIt.lat5,
575  metaPrimitiveIt.lat6,
576  metaPrimitiveIt.lat7,
577  metaPrimitiveIt.lat8};
578 
579  // phiTP (extended DF)
580  outExtP2Ph.emplace_back(
581  L1Phase2MuDTExtPhDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
582  chId.wheel(), // uwh (m_wheel)
583  sectorTP, // usc (m_sector)
584  chId.station(), // ust (m_station)
585  sl, // ust (m_station)
586  (int)round(metaPrimitiveIt.phi * PHIRES_CONV), // uphi (m_phiAngle)
587  (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV), // uphib (m_phiBending)
588  metaPrimitiveIt.quality, // uqua (m_qualityCode)
589  metaPrimitiveIt.index, // uind (m_segmentIndex)
590  (int)round(metaPrimitiveIt.t0) - shift_back * LHC_CLK_FREQ, // ut0 (m_t0Segment)
591  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
592  (int)round(metaPrimitiveIt.x * 1000), // ux (m_xLocal)
593  (int)round(metaPrimitiveIt.tanPhi * 1000), // utan (m_tanPsi)
594  (int)round(metaPrimitiveIt.phi_cmssw * PHIRES_CONV), // uphi (m_phiAngleCMSSW)
595  (int)round(metaPrimitiveIt.phiB_cmssw * PHIBRES_CONV), // uphib (m_phiBendingCMSSW)
596  metaPrimitiveIt.rpcFlag, // urpc (m_rpcFlag)
597  pathWireId,
598  pathTDC,
599  pathLat));
600  }
601  if (df_extended_ == 0 || df_extended_ == 2) {
602  // phiTP (standard DF)
603  outP2Ph.push_back(L1Phase2MuDTPhDigi(
604  (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
605  chId.wheel(), // uwh (m_wheel)
606  sectorTP, // usc (m_sector)
607  chId.station(), // ust (m_station)
608  sl, // ust (m_station)
609  (int)round(metaPrimitiveIt.phi * PHIRES_CONV), // uphi (_phiAngle)
610  (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV), // uphib (m_phiBending)
611  metaPrimitiveIt.quality, // uqua (m_qualityCode)
612  metaPrimitiveIt.index, // uind (m_segmentIndex)
613  (int)round(metaPrimitiveIt.t0) - shift_back * LHC_CLK_FREQ, // ut0 (m_t0Segment)
614  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
615  metaPrimitiveIt.rpcFlag // urpc (m_rpcFlag)
616  ));
617  }
618  } else {
619  if (df_extended_ == 1 || df_extended_ == 2) {
620  int pathWireId[4] = {metaPrimitiveIt.wi1, metaPrimitiveIt.wi2, metaPrimitiveIt.wi3, metaPrimitiveIt.wi4};
621 
622  int pathTDC[4] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
623  max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
624  max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
625  max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1)};
626 
627  int pathLat[4] = {metaPrimitiveIt.lat1, metaPrimitiveIt.lat2, metaPrimitiveIt.lat3, metaPrimitiveIt.lat4};
628 
629  // thTP (extended DF)
630  outExtP2Th.emplace_back(
631  L1Phase2MuDTExtThDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
632  chId.wheel(), // uwh (m_wheel)
633  sectorTP, // usc (m_sector)
634  chId.station(), // ust (m_station)
635  (int)round(metaPrimitiveIt.phi * ZRES_CONV), // uz (m_zGlobal)
636  (int)round(metaPrimitiveIt.phiB * KRES_CONV), // uk (m_kSlope)
637  metaPrimitiveIt.quality, // uqua (m_qualityCode)
638  metaPrimitiveIt.index, // uind (m_segmentIndex)
639  (int)round(metaPrimitiveIt.t0) - shift_back * LHC_CLK_FREQ, // ut0 (m_t0Segment)
640  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
641  (int)round(metaPrimitiveIt.x * 1000), // ux (m_yLocal)
642  (int)round(metaPrimitiveIt.phi_cmssw * PHIRES_CONV), // uphi (m_zCMSSW)
643  (int)round(metaPrimitiveIt.phiB_cmssw * PHIBRES_CONV), // uphib (m_kCMSSW)
644  metaPrimitiveIt.rpcFlag, // urpc (m_rpcFlag)
645  pathWireId,
646  pathTDC,
647  pathLat));
648  }
649  if (df_extended_ == 0 || df_extended_ == 2) {
650  // thTP (standard DF)
651  outP2Th.push_back(L1Phase2MuDTThDigi(
652  (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
653  chId.wheel(), // uwh (m_wheel)
654  sectorTP, // usc (m_sector)
655  chId.station(), // ust (m_station)
656  (int)round(metaPrimitiveIt.phi * ZRES_CONV), // uz (m_zGlobal)
657  (int)round(metaPrimitiveIt.phiB * KRES_CONV), // uk (m_kSlope)
658  metaPrimitiveIt.quality, // uqua (m_qualityCode)
659  metaPrimitiveIt.index, // uind (m_segmentIndex)
660  (int)round(metaPrimitiveIt.t0) - shift_back * LHC_CLK_FREQ, // ut0 (m_t0Segment)
661  (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV), // uchi2 (m_chi2Segment)
662  metaPrimitiveIt.rpcFlag // urpc (m_rpcFlag)
663  ));
664  }
665  }
666  }
667 
668  // Storing RPC hits that were not used elsewhere
669  if (useRPC_) {
670  for (auto rpc_dt_digi = rpc_integrator_->rpcRecHits_translated_.begin();
671  rpc_dt_digi != rpc_integrator_->rpcRecHits_translated_.end();
672  rpc_dt_digi++) {
673  outP2Ph.push_back(*rpc_dt_digi);
674  }
675  }
676 
677  // Storing Phi results
678  if (df_extended_ == 1 || df_extended_ == 2) {
679  std::unique_ptr<L1Phase2MuDTExtPhContainer> resultExtP2Ph(new L1Phase2MuDTExtPhContainer);
680  resultExtP2Ph->setContainer(outExtP2Ph);
681  iEvent.put(std::move(resultExtP2Ph));
682  }
683  if (df_extended_ == 0 || df_extended_ == 2) {
684  std::unique_ptr<L1Phase2MuDTPhContainer> resultP2Ph(new L1Phase2MuDTPhContainer);
685  resultP2Ph->setContainer(outP2Ph);
686  iEvent.put(std::move(resultP2Ph));
687  }
688  outExtP2Ph.clear();
689  outExtP2Ph.erase(outExtP2Ph.begin(), outExtP2Ph.end());
690  outP2Ph.clear();
691  outP2Ph.erase(outP2Ph.begin(), outP2Ph.end());
692 
693  // Storing Theta results
694  if (df_extended_ == 1 || df_extended_ == 2) {
695  std::unique_ptr<L1Phase2MuDTExtThContainer> resultExtP2Th(new L1Phase2MuDTExtThContainer);
696  resultExtP2Th->setContainer(outExtP2Th);
697  iEvent.put(std::move(resultExtP2Th));
698  }
699  if (df_extended_ == 0 || df_extended_ == 2) {
700  std::unique_ptr<L1Phase2MuDTThContainer> resultP2Th(new L1Phase2MuDTThContainer);
701  resultP2Th->setContainer(outP2Th);
702  iEvent.put(std::move(resultP2Th));
703  }
704  outExtP2Th.clear();
705  outExtP2Th.erase(outExtP2Th.begin(), outExtP2Th.end());
706  outP2Th.clear();
707  outP2Th.erase(outP2Th.begin(), outP2Th.end());
708 }
709 
710 void DTTrigPhase2Prod::endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
711  grouping_obj_->finish();
712  mpathanalyzer_->finish();
713  mpathqualityenhancer_->finish();
714  mpathredundantfilter_->finish();
715  mpathassociator_->finish();
716  rpc_integrator_->finish();
717 };
718 
719 bool DTTrigPhase2Prod::outer(const metaPrimitive& mp) const {
720  int counter = (mp.wi5 != -1) + (mp.wi6 != -1) + (mp.wi7 != -1) + (mp.wi8 != -1);
721  return (counter > 2);
722 }
723 
724 bool DTTrigPhase2Prod::inner(const metaPrimitive& mp) const {
725  int counter = (mp.wi1 != -1) + (mp.wi2 != -1) + (mp.wi3 != -1) + (mp.wi4 != -1);
726  return (counter > 2);
727 }
728 
729 bool DTTrigPhase2Prod::hasPosRF(int wh, int sec) const { return wh > 0 || (wh == 0 && sec % 4 > 1); }
730 
731 void DTTrigPhase2Prod::printmP(const string& ss, const metaPrimitive& mP) const {
732  DTSuperLayerId slId(mP.rawId);
733  LogInfo("DTTrigPhase2Prod") << ss << (int)slId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
734  << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
735  << setw(5) << left << mP.tdc1 << " " << setw(5) << left << mP.tdc2 << " " << setw(5)
736  << left << mP.tdc3 << " " << setw(5) << left << mP.tdc4 << " " << setw(10) << right
737  << mP.x << " " << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " "
738  << setw(13) << left << mP.chi2 << " r:" << rango(mP);
739 }
740 
741 void DTTrigPhase2Prod::printmPC(const string& ss, const metaPrimitive& mP) const {
742  DTChamberId ChId(mP.rawId);
743  LogInfo("DTTrigPhase2Prod") << ss << (int)ChId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
744  << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
745  << setw(2) << left << mP.wi5 << " " << setw(2) << left << mP.wi6 << " " << setw(2) << left
746  << mP.wi7 << " " << setw(2) << left << mP.wi8 << " " << setw(5) << left << mP.tdc1 << " "
747  << setw(5) << left << mP.tdc2 << " " << setw(5) << left << mP.tdc3 << " " << setw(5)
748  << left << mP.tdc4 << " " << setw(5) << left << mP.tdc5 << " " << setw(5) << left
749  << mP.tdc6 << " " << setw(5) << left << mP.tdc7 << " " << setw(5) << left << mP.tdc8
750  << " " << setw(2) << left << mP.lat1 << " " << setw(2) << left << mP.lat2 << " "
751  << setw(2) << left << mP.lat3 << " " << setw(2) << left << mP.lat4 << " " << setw(2)
752  << left << mP.lat5 << " " << setw(2) << left << mP.lat6 << " " << setw(2) << left
753  << mP.lat7 << " " << setw(2) << left << mP.lat8 << " " << setw(10) << right << mP.x << " "
754  << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13)
755  << left << mP.chi2 << " r:" << rango(mP);
756 }
757 
759  if (mp.quality == 1 or mp.quality == 2)
760  return 3;
761  if (mp.quality == 3 or mp.quality == 4)
762  return 4;
763  return mp.quality;
764 }
765 
766 void DTTrigPhase2Prod::assignIndex(std::vector<metaPrimitive>& inMPaths) {
767  std::map<int, std::vector<metaPrimitive>> primsPerBX;
768  for (const auto& metaPrimitive : inMPaths) {
769  int BX = round(metaPrimitive.t0 / 25.);
770  primsPerBX[BX].push_back(metaPrimitive);
771  }
772  inMPaths.clear();
773  for (auto& prims : primsPerBX) {
774  assignIndexPerBX(prims.second);
775  for (const auto& primitive : prims.second)
776  inMPaths.push_back(primitive);
777  }
778 }
779 
780 void DTTrigPhase2Prod::assignIndexPerBX(std::vector<metaPrimitive>& inMPaths) {
781  // First we asociate a new index to the metaprimitive depending on quality or phiB;
782  uint32_t rawId = -1;
783  int numP = -1;
784  for (auto& metaPrimitiveIt : inMPaths) {
785  numP++;
786  rawId = metaPrimitiveIt.rawId;
787  int iOrder = assignQualityOrder(metaPrimitiveIt);
788  int inf = 0;
789  int numP2 = -1;
790  for (auto& metaPrimitiveItN : inMPaths) {
791  int nOrder = assignQualityOrder(metaPrimitiveItN);
792  numP2++;
793  if (rawId != metaPrimitiveItN.rawId)
794  continue;
795  if (numP2 == numP) {
796  metaPrimitiveIt.index = inf;
797  break;
798  } else if (iOrder < nOrder) {
799  inf++;
800  } else if (iOrder > nOrder) {
801  metaPrimitiveItN.index++;
802  } else if (iOrder == nOrder) {
803  if (std::abs(metaPrimitiveIt.phiB) >= std::abs(metaPrimitiveItN.phiB)) {
804  inf++;
805  } else if (std::abs(metaPrimitiveIt.phiB) < std::abs(metaPrimitiveItN.phiB)) {
806  metaPrimitiveItN.index++;
807  }
808  }
809  } // ending second for
810  } // ending first for
811 }
812 
814  if (mP.quality > 9 || mP.quality < 1)
815  return -1;
816 
817  return qmap_.find(mP.quality)->second;
818 }
819 
820 std::vector<DTDigiCollection*> DTTrigPhase2Prod::distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ) {
821  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*> tmpVector;
822  tmpVector.clear();
823  std::vector<DTDigiCollection*> collVector;
824  collVector.clear();
825  while (!inQ.empty()) {
826  processDigi(inQ, tmpVector);
827  }
828  for (auto& sQ : tmpVector) {
829  DTDigiCollection tmpColl;
830  while (!sQ->empty()) {
831  tmpColl.insertDigi((sQ->front().first), (sQ->front().second));
832  sQ->pop();
833  }
834  collVector.push_back(&tmpColl);
835  }
836  return collVector;
837 }
838 
839 void DTTrigPhase2Prod::processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
840  std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec) {
841  bool classified = false;
842  if (!vec.empty()) {
843  for (auto& sC : vec) { // Conditions for entering a super cell.
844  if ((sC->front().second.time() + superCelltimewidth_) > inQ.front().second.time()) {
845  // Time requirement
846  if (TMath::Abs(sC->front().second.wire() - inQ.front().second.wire()) <= superCellhalfspacewidth_) {
847  // Spatial requirement
848  sC->push(std::move(inQ.front()));
849  classified = true;
850  }
851  }
852  }
853  }
854  if (classified) {
855  inQ.pop();
856  return;
857  }
858 
859  std::queue<std::pair<DTLayerId, DTDigi>> newQueue;
860 
861  std::pair<DTLayerId, DTDigi> tmpPair;
862  tmpPair = std::move(inQ.front());
863  newQueue.push(tmpPair);
864  inQ.pop();
865 
866  vec.push_back(&newQueue);
867 }
868 
int assignQualityOrder(const metaPrimitive &mP) const
EventNumber_t event() const
Definition: EventID.h:40
std::unique_ptr< MotherGrouping > grouping_obj_
std::map< DTChamberId, DTDigiCollection, std::less< DTChamberId > > DTDigiMap
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool hasPosRF(int wh, int sec) const
RunID const & id() const
Definition: RunBase.h:39
void produce(edm::Event &iEvent, const edm::EventSetup &iEventSetup) override
Producer: process every event and generates trigger data.
edm::ConsumesCollector consumesColl(consumesCollector())
constexpr int CHI2RES_CONV
Definition: constants.h:226
constexpr float PHIBRES_CONV
Definition: constants.h:225
RunNumber_t run() const
Definition: RunID.h:36
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
mpathqualityenhancer_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
mpathassociator_
DTTrigPhase2Prod(const edm::ParameterSet &pset)
Constructor.
DTChamberId chamberId() const
Return the corresponding ChamberId.
dtDigisToken_
DTGeometry const * dtGeo_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
superCellhalfspacewidth_
std::vector< MuonPathPtr > MuonPathPtrs
Definition: MuonPath.h:128
geometry_tag_
produces< L1Phase2MuDTExtThContainer >()
edm::EDGetTokenT< DTDigiCollection > dtDigisToken_
void processDigi(std::queue< std::pair< DTLayerId, DTDigi >> &inQ, std::vector< std::queue< std::pair< DTLayerId, DTDigi >> * > &vec)
std::vector< DTDigiCollection * > distribDigis(std::queue< std::pair< DTLayerId, DTDigi >> &inQ)
std::unique_ptr< RPCIntegrator > rpc_integrator_
void printmPC(const std::string &ss, const metaPrimitive &mP) const
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:45
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomH
mpathanalyzer_
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsLabel_
Primitive< F, X >::type primitive(const F &f)
Definition: Primitive.h:41
const uint16_t range(const Frame &aFrame)
produces< L1Phase2MuDTExtPhContainer >()
void assignIndex(std::vector< metaPrimitive > &inMPaths)
int iEvent
Definition: GenABIO.cc:224
globalcoordsobtainer_
do_correlation_
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:32
std::unique_ptr< MuonPathAnalyzer > mpathanalyzer_
constexpr float KRES_CONV
Definition: constants.h:229
def move
Definition: eostools.py:511
string inf
Definition: EcalCondDB.py:96
std::unique_ptr< MPFilter > mpathqualityenhancer_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int superLayer() const
Return the superlayer number.
std::unique_ptr< MuonPathAssociator > mpathassociator_
constexpr float ZRES_CONV
Definition: constants.h:228
bool inner(const metaPrimitive &mp) const
rpc_integrator_
mpathredundantfilter_
Log< level::Info, false > LogInfo
MP_QUALITY
Definition: constants.h:42
void endRun(edm::Run const &iRun, const edm::EventSetup &iEventSetup) override
endRun: finish things
std::unique_ptr< MPFilter > mpathredundantfilter_
std::shared_ptr< GlobalCoordsObtainer > globalcoordsobtainer_
std::pair< const_iterator, const_iterator > Range
bool outer(const metaPrimitive &mp) const
std::string geometry_tag_
double b
Definition: hdecay.h:118
std::vector< DigiType >::const_iterator const_iterator
constexpr int LHC_CLK_FREQ
Definition: constants.h:176
DTDigiMap::iterator DTDigiMap_iterator
superCelltimewidth_
edm::EventID id() const
Definition: EventBase.h:59
void assignIndexPerBX(std::vector< metaPrimitive > &inMPaths)
~DTTrigPhase2Prod() override
Destructor.
int sector() const
Definition: DTChamberId.h:49
double a
Definition: hdecay.h:119
static std::atomic< unsigned int > counter
T get() const
Definition: EventSetup.h:82
DTDigiMap::const_iterator DTDigiMap_const_iterator
constexpr float PHIRES_CONV
Definition: constants.h:224
activateBuffer_
df_extended_
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
const std::unordered_map< int, int > qmap_
void printmP(const std::string &ss, const metaPrimitive &mP) const
void beginRun(edm::Run const &iRun, const edm::EventSetup &iEventSetup) override
Create Trigger Units before starting event processing.
int rango(const metaPrimitive &mp) const
std::vector< std::pair< int, MuonPath > > primitives_
rpcRecHitsLabel_
Definition: Run.h:45
produces< L1Phase2MuDTThContainer >()
#define LogDebug(id)