Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:46

0001 #include "FWCore/Framework/interface/MakerMacros.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/Framework/interface/ESHandle.h"
0004 #include "FWCore/Utilities/interface/ESGetToken.h"
0005 #include "FWCore/Framework/interface/ModuleFactory.h"
0006 #include "FWCore/Framework/interface/ESProducer.h"
0007 #include "FWCore/Framework/interface/ESHandle.h"
0008 #include "FWCore/Framework/interface/ESProducts.h"
0009 #include "FWCore/Framework/interface/ConsumesCollector.h"
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/Framework/interface/Run.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 
0018 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0020 #include "Geometry/DTGeometry/interface/DTLayer.h"
0021 
0022 #include "L1Trigger/DTTriggerPhase2/interface/MuonPath.h"
0023 #include "L1Trigger/DTTriggerPhase2/interface/constants.h"
0024 
0025 #include "L1Trigger/DTTriggerPhase2/interface/MotherGrouping.h"
0026 #include "L1Trigger/DTTriggerPhase2/interface/TrapezoidalGrouping.h"
0027 #include "L1Trigger/DTTriggerPhase2/interface/HoughGrouping.h"
0028 #include "L1Trigger/DTTriggerPhase2/interface/PseudoBayesGrouping.h"
0029 #include "L1Trigger/DTTriggerPhase2/interface/LateralityProvider.h"
0030 #include "L1Trigger/DTTriggerPhase2/interface/LateralityBasicProvider.h"
0031 #include "L1Trigger/DTTriggerPhase2/interface/LateralityCoarsedProvider.h"
0032 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathAnalyzer.h"
0033 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathSLFitter.h"
0034 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathCorFitter.h"
0035 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathAnalyticAnalyzer.h"
0036 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathAnalyzerInChamber.h"
0037 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathAssociator.h"
0038 #include "L1Trigger/DTTriggerPhase2/interface/MuonPathConfirmator.h"
0039 #include "L1Trigger/DTTriggerPhase2/interface/MPFilter.h"
0040 #include "L1Trigger/DTTriggerPhase2/interface/MPSLFilter.h"
0041 #include "L1Trigger/DTTriggerPhase2/interface/MPCorFilter.h"
0042 #include "L1Trigger/DTTriggerPhase2/interface/MPQualityEnhancerFilter.h"
0043 #include "L1Trigger/DTTriggerPhase2/interface/MPRedundantFilter.h"
0044 #include "L1Trigger/DTTriggerPhase2/interface/MPCleanHitsFilter.h"
0045 #include "L1Trigger/DTTriggerPhase2/interface/MPQualityEnhancerFilterBayes.h"
0046 #include "L1Trigger/DTTriggerPhase2/interface/GlobalCoordsObtainer.h"
0047 
0048 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0049 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
0050 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0051 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0052 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
0053 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTPhContainer.h"
0054 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTPhDigi.h"
0055 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtPhContainer.h"
0056 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtPhDigi.h"
0057 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTThContainer.h"
0058 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTThDigi.h"
0059 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtThContainer.h"
0060 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtThDigi.h"
0061 
0062 // DT trigger GeomUtils
0063 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
0064 
0065 //RPC TP
0066 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
0067 #include <DataFormats/MuonDetId/interface/RPCDetId.h>
0068 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0069 #include "L1Trigger/DTTriggerPhase2/interface/RPCIntegrator.h"
0070 
0071 #include <fstream>
0072 #include <iostream>
0073 #include <queue>
0074 #include <cmath>
0075 
0076 using namespace edm;
0077 using namespace std;
0078 using namespace cmsdt;
0079 
0080 class DTTrigPhase2Prod : public edm::stream::EDProducer<> {
0081   typedef std::map<DTChamberId, DTDigiCollection, std::less<DTChamberId>> DTDigiMap;
0082   typedef DTDigiMap::iterator DTDigiMap_iterator;
0083   typedef DTDigiMap::const_iterator DTDigiMap_const_iterator;
0084 
0085 public:
0086   //! Constructor
0087   DTTrigPhase2Prod(const edm::ParameterSet& pset);
0088 
0089   //! Destructor
0090   ~DTTrigPhase2Prod() override;
0091 
0092   //! Create Trigger Units before starting event processing
0093   void beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
0094 
0095   //! Producer: process every event and generates trigger data
0096   void produce(edm::Event& iEvent, const edm::EventSetup& iEventSetup) override;
0097 
0098   //! endRun: finish things
0099   void endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
0100 
0101   // Methods
0102   int rango(const metaPrimitive& mp) const;
0103   bool outer(const metaPrimitive& mp) const;
0104   bool inner(const metaPrimitive& mp) const;
0105   void printmP(const std::string& ss, const metaPrimitive& mP) const;
0106   void printmP(const metaPrimitive& mP) const;
0107   void printmPC(const std::string& ss, const metaPrimitive& mP) const;
0108   void printmPC(const metaPrimitive& mP) const;
0109   bool hasPosRF(int wh, int sec) const;
0110 
0111   // Getter-methods
0112   MP_QUALITY getMinimumQuality(void);
0113 
0114   // Setter-methods
0115   void setChiSquareThreshold(float ch2Thr);
0116   void setMinimumQuality(MP_QUALITY q);
0117 
0118   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0119 
0120   // data-members
0121   const DTGeometry* dtGeo_;
0122   edm::ESGetToken<DTGeometry, MuonGeometryRecord> dtGeomH;
0123   std::vector<std::pair<int, MuonPath>> primitives_;
0124 
0125 private:
0126   // Trigger Configuration Manager CCB validity flag
0127   bool my_CCBValid_;
0128 
0129   // BX offset used to correct DTTPG output
0130   int my_BXoffset_;
0131 
0132   // Debug Flag
0133   bool debug_;
0134   bool dump_;
0135   double dT0_correlate_TP_;
0136   int scenario_;
0137   int df_extended_;
0138   int max_index_;
0139 
0140   bool output_mixer_;
0141   bool output_latpredictor_;
0142   bool output_slfitter_;
0143   bool output_slfilter_;
0144   bool output_confirmed_;
0145   bool output_matcher_;
0146   bool skip_processing_;
0147   bool allow_confirmation_;
0148 
0149   // ParameterSet
0150   edm::EDGetTokenT<DTDigiCollection> dtDigisToken_;
0151   edm::EDGetTokenT<RPCRecHitCollection> rpcRecHitsLabel_;
0152 
0153   // Grouping attributes and methods
0154   int algo_;  // Grouping code
0155   std::unique_ptr<MotherGrouping> grouping_obj_;
0156   std::unique_ptr<MuonPathAnalyzer> mpathanalyzer_;
0157   std::unique_ptr<LateralityProvider> latprovider_;
0158   std::unique_ptr<MPFilter> mpathqualityenhancer_;
0159   std::unique_ptr<MPFilter> mpathqualityenhancerbayes_;
0160   std::unique_ptr<MPFilter> mpathredundantfilter_;
0161   std::unique_ptr<MPFilter> mpathhitsfilter_;
0162   std::unique_ptr<MuonPathAnalyzer> mpathassociator_;
0163   std::unique_ptr<MuonPathConfirmator> mpathconfirmator_;
0164   std::unique_ptr<MPFilter> mpathcorfilter_;
0165   std::shared_ptr<GlobalCoordsObtainer> globalcoordsobtainer_;
0166 
0167   // Buffering
0168   bool activateBuffer_;
0169   int superCellhalfspacewidth_;
0170   float superCelltimewidth_;
0171   std::vector<DTDigiCollection*> distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ);
0172   void processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
0173                    std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec);
0174 
0175   // RPC
0176   std::unique_ptr<RPCIntegrator> rpc_integrator_;
0177   bool useRPC_;
0178 
0179   void assignIndex(std::vector<metaPrimitive>& inMPaths);
0180   void assignIndexPerBX(std::vector<metaPrimitive>& inMPaths);
0181   int assignQualityOrder(const metaPrimitive& mP) const;
0182 
0183   const std::unordered_map<int, int> qmap_;
0184 };
0185 
0186 namespace {
0187   struct {
0188     bool operator()(std::pair<DTLayerId, DTDigi> a, std::pair<DTLayerId, DTDigi> b) const {
0189       return (a.second.time() < b.second.time());
0190     }
0191   } const DigiTimeOrdering;
0192 }  // namespace
0193 
0194 DTTrigPhase2Prod::DTTrigPhase2Prod(const ParameterSet& pset)
0195     : qmap_({{8, 8}, {7, 7}, {6, 6}, {4, 4}, {3, 3}, {2, 2}, {1, 1}}) {
0196   produces<L1Phase2MuDTPhContainer>();
0197   produces<L1Phase2MuDTThContainer>();
0198   produces<L1Phase2MuDTExtPhContainer>();
0199   produces<L1Phase2MuDTExtThContainer>();
0200 
0201   debug_ = pset.getUntrackedParameter<bool>("debug");
0202   dump_ = pset.getUntrackedParameter<bool>("dump");
0203 
0204   scenario_ = pset.getParameter<int>("scenario");
0205 
0206   df_extended_ = pset.getParameter<int>("df_extended");
0207   max_index_ = pset.getParameter<int>("max_primitives") - 1;
0208 
0209   dtDigisToken_ = consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiTag"));
0210 
0211   rpcRecHitsLabel_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("rpcRecHits"));
0212   useRPC_ = pset.getParameter<bool>("useRPC");
0213 
0214   // Choosing grouping scheme:
0215   algo_ = pset.getParameter<int>("algo");
0216 
0217   // shortcuts
0218 
0219   output_mixer_ = pset.getParameter<bool>("output_mixer");
0220   output_latpredictor_ = pset.getParameter<bool>("output_latpredictor");
0221   output_slfitter_ = pset.getParameter<bool>("output_slfitter");
0222   output_slfilter_ = pset.getParameter<bool>("output_slfilter");
0223   output_confirmed_ = pset.getParameter<bool>("output_confirmed");
0224   output_matcher_ = pset.getParameter<bool>("output_matcher");
0225   allow_confirmation_ = pset.getParameter<bool>("allow_confirmation");
0226 
0227   edm::ConsumesCollector consumesColl(consumesCollector());
0228   globalcoordsobtainer_ = std::make_shared<GlobalCoordsObtainer>(pset);
0229   globalcoordsobtainer_->generate_luts();
0230 
0231   if (algo_ == PseudoBayes) {
0232     grouping_obj_ =
0233         std::make_unique<PseudoBayesGrouping>(pset.getParameter<edm::ParameterSet>("PseudoBayesPattern"), consumesColl);
0234   } else if (algo_ == HoughTrans) {
0235     grouping_obj_ =
0236         std::make_unique<HoughGrouping>(pset.getParameter<edm::ParameterSet>("HoughGrouping"), consumesColl);
0237   } else {
0238     grouping_obj_ = std::make_unique<TrapezoidalGrouping>(pset, consumesColl);
0239   }
0240 
0241   if (algo_ == Standard) {
0242     if (debug_)
0243       LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: JM analyzer";
0244     mpathanalyzer_ = std::make_unique<MuonPathSLFitter>(pset, consumesColl, globalcoordsobtainer_);
0245     latprovider_ = std::make_unique<LateralityCoarsedProvider>(pset, consumesColl);
0246   } else {
0247     if (debug_)
0248       LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: Full chamber analyzer";
0249     mpathanalyzer_ = std::make_unique<MuonPathAnalyzerInChamber>(pset, consumesColl, globalcoordsobtainer_);
0250   }
0251 
0252   // Getting buffer option
0253   activateBuffer_ = pset.getParameter<bool>("activateBuffer");
0254   superCellhalfspacewidth_ = pset.getParameter<int>("superCellspacewidth") / 2;
0255   superCelltimewidth_ = pset.getParameter<double>("superCelltimewidth");
0256 
0257   mpathqualityenhancer_ = std::make_unique<MPSLFilter>(pset);
0258   mpathqualityenhancerbayes_ = std::make_unique<MPQualityEnhancerFilterBayes>(pset);
0259   mpathredundantfilter_ = std::make_unique<MPRedundantFilter>(pset);
0260   mpathhitsfilter_ = std::make_unique<MPCleanHitsFilter>(pset);
0261   mpathconfirmator_ = std::make_unique<MuonPathConfirmator>(pset, consumesColl);
0262   mpathassociator_ = std::make_unique<MuonPathCorFitter>(pset, consumesColl, globalcoordsobtainer_);
0263   mpathcorfilter_ = std::make_unique<MPCorFilter>(pset);
0264   rpc_integrator_ = std::make_unique<RPCIntegrator>(pset, consumesColl);
0265 
0266   dtGeomH = esConsumes<DTGeometry, MuonGeometryRecord, edm::Transition::BeginRun>();
0267 }
0268 
0269 DTTrigPhase2Prod::~DTTrigPhase2Prod() {
0270   if (debug_)
0271     LogDebug("DTTrigPhase2Prod") << "DTp2: calling destructor" << std::endl;
0272 }
0273 
0274 void DTTrigPhase2Prod::beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
0275   if (debug_)
0276     LogDebug("DTTrigPhase2Prod") << "beginRun " << iRun.id().run();
0277   if (debug_)
0278     LogDebug("DTTrigPhase2Prod") << "beginRun: getting DT geometry";
0279 
0280   grouping_obj_->initialise(iEventSetup);               // Grouping object initialisation
0281   mpathanalyzer_->initialise(iEventSetup);              // Analyzer object initialisation
0282   mpathqualityenhancer_->initialise(iEventSetup);       // Filter object initialisation
0283   mpathredundantfilter_->initialise(iEventSetup);       // Filter object initialisation
0284   mpathqualityenhancerbayes_->initialise(iEventSetup);  // Filter object initialisation
0285   mpathhitsfilter_->initialise(iEventSetup);
0286   mpathassociator_->initialise(iEventSetup);  // Associator object initialisation
0287   mpathcorfilter_->initialise(iEventSetup);
0288 
0289   if (auto geom = iEventSetup.getHandle(dtGeomH)) {
0290     dtGeo_ = &(*geom);
0291   }
0292 }
0293 
0294 void DTTrigPhase2Prod::produce(Event& iEvent, const EventSetup& iEventSetup) {
0295   if (debug_)
0296     LogDebug("DTTrigPhase2Prod") << "produce";
0297   edm::Handle<DTDigiCollection> dtdigis;
0298   iEvent.getByToken(dtDigisToken_, dtdigis);
0299 
0300   if (debug_)
0301     LogDebug("DTTrigPhase2Prod") << "\t Getting the RPC RecHits" << std::endl;
0302   edm::Handle<RPCRecHitCollection> rpcRecHits;
0303   iEvent.getByToken(rpcRecHitsLabel_, rpcRecHits);
0304 
0305   ////////////////////////////////
0306   // GROUPING CODE:
0307   ////////////////////////////////
0308 
0309   DTDigiMap digiMap;
0310   DTDigiCollection::DigiRangeIterator detUnitIt;
0311   for (const auto& detUnitIt : *dtdigis) {
0312     const DTLayerId& layId = detUnitIt.first;
0313     const DTChamberId chambId = layId.superlayerId().chamberId();
0314     const DTDigiCollection::Range& range = detUnitIt.second;
0315     digiMap[chambId].put(range, layId);
0316   }
0317 
0318   // generate a list muon paths for each event!!!
0319   if (debug_ && activateBuffer_)
0320     LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber using a buffer and super cells.";
0321   else if (debug_)
0322     LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber.";
0323 
0324   std::map<int, MuonPathPtrs> muonpaths;
0325   for (const auto& ich : dtGeo_->chambers()) {
0326     // 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.
0327     const DTChamber* chamb = ich;
0328     DTChamberId chid = chamb->id();
0329     DTDigiMap_iterator dmit = digiMap.find(chid);
0330 
0331     if (dmit == digiMap.end())
0332       continue;
0333 
0334     if (activateBuffer_) {  // Use buffering (per chamber) or not
0335       // Import digis from the station
0336       std::vector<std::pair<DTLayerId, DTDigi>> tmpvec;
0337       tmpvec.clear();
0338 
0339       for (const auto& dtLayerIdIt : (*dmit).second) {
0340         for (DTDigiCollection::const_iterator digiIt = (dtLayerIdIt.second).first;
0341              digiIt != (dtLayerIdIt.second).second;
0342              digiIt++) {
0343           tmpvec.emplace_back(dtLayerIdIt.first, *digiIt);
0344         }
0345       }
0346 
0347       // Check to enhance CPU time usage
0348       if (tmpvec.empty())
0349         continue;
0350 
0351       // Order digis depending on TDC time and insert them into a queue (FIFO buffer). TODO: adapt for MC simulations.
0352       std::sort(tmpvec.begin(), tmpvec.end(), DigiTimeOrdering);
0353       std::queue<std::pair<DTLayerId, DTDigi>> timequeue;
0354 
0355       for (const auto& elem : tmpvec)
0356         timequeue.emplace(elem);
0357       tmpvec.clear();
0358 
0359       // Distribute the digis from the queue into supercells
0360       std::vector<DTDigiCollection*> superCells;
0361       superCells = distribDigis(timequeue);
0362 
0363       // Process each supercell & collect the resulting muonpaths (as the muonpaths std::vector is only enlarged each time
0364       // the groupings access it, it's not needed to "collect" the final products).
0365 
0366       while (!superCells.empty()) {
0367         grouping_obj_->run(iEvent, iEventSetup, *(superCells.back()), muonpaths[chid.rawId()]);
0368         superCells.pop_back();
0369       }
0370     } else {
0371       grouping_obj_->run(iEvent, iEventSetup, (*dmit).second, muonpaths[chid.rawId()]);
0372     }
0373   }
0374   digiMap.clear();
0375 
0376   if (dump_) {
0377     for (auto& ch_muonpaths : muonpaths) {
0378       for (unsigned int i = 0; i < ch_muonpaths.second.size(); i++) {
0379         stringstream ss;
0380         ss << iEvent.id().event() << "      mpath " << i << ": ";
0381         for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0382           ss << ch_muonpaths.second.at(i)->primitive(lay)->channelId() << " ";
0383         for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0384           ss << ch_muonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
0385         for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0386           ss << ch_muonpaths.second.at(i)->primitive(lay)->laterality() << " ";
0387         LogInfo("DTTrigPhase2Prod") << ss.str();
0388       }
0389     }
0390   }
0391 
0392   std::map<int, std::vector<lat_vector>> lateralities;
0393   if (!output_mixer_) {
0394     for (auto& ch_muonpaths : muonpaths) {
0395       if (algo_ == Standard) {
0396         latprovider_->run(iEvent, iEventSetup, ch_muonpaths.second, lateralities[ch_muonpaths.first]);
0397       }
0398     }
0399   }
0400 
0401   // FILTER GROUPING
0402   std::map<int, MuonPathPtrs> filteredmuonpaths;
0403   for (auto& ch_muonpaths : muonpaths) {
0404     if (algo_ == Standard) {
0405       mpathredundantfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
0406     } else {
0407       mpathhitsfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
0408     }
0409   }
0410 
0411   if (dump_) {
0412     for (auto& ch_filteredmuonpaths : filteredmuonpaths) {
0413       for (unsigned int i = 0; i < ch_filteredmuonpaths.second.size(); i++) {
0414         stringstream ss;
0415         ss << iEvent.id().event() << " filt. mpath " << i << ": ";
0416         for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
0417           ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->channelId() << " ";
0418         for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
0419           ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
0420         LogInfo("DTTrigPhase2Prod") << ss.str();
0421       }
0422     }
0423   }
0424 
0425   skip_processing_ = output_mixer_ || output_latpredictor_;
0426 
0427   ///////////////////////////////////////////
0428   /// Fitting SECTION;
0429   ///////////////////////////////////////////
0430 
0431   if (debug_) {
0432     for (auto& ch_muonpaths : muonpaths) {
0433       LogDebug("DTTrigPhase2Prod") << "MUON PATHS found: " << ch_muonpaths.second.size() << " ("
0434                                    << filteredmuonpaths[ch_muonpaths.first].size() << ") in event "
0435                                    << iEvent.id().event();
0436     }
0437   }
0438   if (debug_)
0439     LogDebug("DTTrigPhase2Prod") << "filling NmetaPrimtives" << std::endl;
0440   std::map<int, std::vector<metaPrimitive>> metaPrimitives;
0441   std::map<int, MuonPathPtrs> outmpaths;
0442   if (algo_ == Standard) {
0443     if (debug_)
0444       LogDebug("DTTrigPhase2Prod") << "Fitting 1SL ";
0445     for (auto& ch_muonpaths : muonpaths) {  // FIXME, do we need filtered muonpaths?
0446       if (!output_mixer_ && !output_latpredictor_)
0447         mpathanalyzer_->run(iEvent,
0448                             iEventSetup,
0449                             ch_muonpaths.second,
0450                             lateralities[ch_muonpaths.first],
0451                             metaPrimitives[ch_muonpaths.first]);
0452       else if (output_mixer_) {
0453         for (auto& inMPath : ch_muonpaths.second) {
0454           auto sl = inMPath->primitive(0)->superLayerId();  // 0, 1, 2
0455           int selected_lay = 1;
0456           if (inMPath->primitive(0)->tdcTimeStamp() != -1)
0457             selected_lay = 0;
0458           int dumLayId = inMPath->primitive(selected_lay)->cameraId();
0459           auto dtDumlayerId = DTLayerId(dumLayId);
0460           DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
0461           if (sl == 0)
0462             metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0463                                                                            -1,
0464                                                                            -1,
0465                                                                            -1,
0466                                                                            -1,
0467                                                                            -1,
0468                                                                            -1,
0469                                                                            -1,
0470                                                                            -1,
0471                                                                            -1,
0472                                                                            inMPath->primitive(0)->channelId(),
0473                                                                            inMPath->primitive(0)->tdcTimeStamp(),
0474                                                                            -1,
0475                                                                            inMPath->primitive(1)->channelId(),
0476                                                                            inMPath->primitive(1)->tdcTimeStamp(),
0477                                                                            -1,
0478                                                                            inMPath->primitive(2)->channelId(),
0479                                                                            inMPath->primitive(2)->tdcTimeStamp(),
0480                                                                            -1,
0481                                                                            inMPath->primitive(3)->channelId(),
0482                                                                            inMPath->primitive(3)->tdcTimeStamp(),
0483                                                                            -1,
0484                                                                            -1,
0485                                                                            -1,
0486                                                                            -1,
0487                                                                            -1,
0488                                                                            -1,
0489                                                                            -1,
0490                                                                            -1,
0491                                                                            -1,
0492                                                                            -1,
0493                                                                            -1,
0494                                                                            -1,
0495                                                                            -1,
0496                                                                            -1}));
0497           else
0498             metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0499                                                                            -1,
0500                                                                            -1,
0501                                                                            -1,
0502                                                                            -1,
0503                                                                            -1,
0504                                                                            -1,
0505                                                                            -1,
0506                                                                            -1,
0507                                                                            -1,
0508                                                                            -1,
0509                                                                            -1,
0510                                                                            -1,
0511                                                                            -1,
0512                                                                            -1,
0513                                                                            -1,
0514                                                                            -1,
0515                                                                            -1,
0516                                                                            -1,
0517                                                                            -1,
0518                                                                            -1,
0519                                                                            -1,
0520                                                                            inMPath->primitive(0)->channelId(),
0521                                                                            inMPath->primitive(0)->tdcTimeStamp(),
0522                                                                            -1,
0523                                                                            inMPath->primitive(1)->channelId(),
0524                                                                            inMPath->primitive(1)->tdcTimeStamp(),
0525                                                                            -1,
0526                                                                            inMPath->primitive(2)->channelId(),
0527                                                                            inMPath->primitive(2)->tdcTimeStamp(),
0528                                                                            -1,
0529                                                                            inMPath->primitive(3)->channelId(),
0530                                                                            inMPath->primitive(3)->tdcTimeStamp(),
0531                                                                            -1,
0532                                                                            -1}));
0533         }
0534       } else if (output_latpredictor_) {
0535         int imp = -1;
0536         for (auto& inMPath : ch_muonpaths.second) {
0537           imp++;
0538           auto sl = inMPath->primitive(0)->superLayerId();  // 0, 1, 2
0539           int selected_lay = 1;
0540           if (inMPath->primitive(0)->tdcTimeStamp() != -1)
0541             selected_lay = 0;
0542           int dumLayId = inMPath->primitive(selected_lay)->cameraId();
0543           auto dtDumlayerId = DTLayerId(dumLayId);
0544           DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
0545           for (auto& latcomb : lateralities[ch_muonpaths.first][imp]) {
0546             if (sl == 0)
0547               metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0548                                                                              -1,
0549                                                                              -1,
0550                                                                              -1,
0551                                                                              -1,
0552                                                                              -1,
0553                                                                              -1,
0554                                                                              -1,
0555                                                                              -1,
0556                                                                              -1,
0557                                                                              inMPath->primitive(0)->channelId(),
0558                                                                              inMPath->primitive(0)->tdcTimeStamp(),
0559                                                                              latcomb[0],
0560                                                                              inMPath->primitive(1)->channelId(),
0561                                                                              inMPath->primitive(1)->tdcTimeStamp(),
0562                                                                              latcomb[1],
0563                                                                              inMPath->primitive(2)->channelId(),
0564                                                                              inMPath->primitive(2)->tdcTimeStamp(),
0565                                                                              latcomb[2],
0566                                                                              inMPath->primitive(3)->channelId(),
0567                                                                              inMPath->primitive(3)->tdcTimeStamp(),
0568                                                                              latcomb[3],
0569                                                                              -1,
0570                                                                              -1,
0571                                                                              -1,
0572                                                                              -1,
0573                                                                              -1,
0574                                                                              -1,
0575                                                                              -1,
0576                                                                              -1,
0577                                                                              -1,
0578                                                                              -1,
0579                                                                              -1,
0580                                                                              -1,
0581                                                                              -1}));
0582             else
0583               metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0584                                                                              -1,
0585                                                                              -1,
0586                                                                              -1,
0587                                                                              -1,
0588                                                                              -1,
0589                                                                              -1,
0590                                                                              -1,
0591                                                                              -1,
0592                                                                              -1,
0593                                                                              -1,
0594                                                                              -1,
0595                                                                              -1,
0596                                                                              -1,
0597                                                                              -1,
0598                                                                              -1,
0599                                                                              -1,
0600                                                                              -1,
0601                                                                              -1,
0602                                                                              -1,
0603                                                                              -1,
0604                                                                              -1,
0605                                                                              inMPath->primitive(0)->channelId(),
0606                                                                              inMPath->primitive(0)->tdcTimeStamp(),
0607                                                                              latcomb[0],
0608                                                                              inMPath->primitive(1)->channelId(),
0609                                                                              inMPath->primitive(1)->tdcTimeStamp(),
0610                                                                              latcomb[1],
0611                                                                              inMPath->primitive(2)->channelId(),
0612                                                                              inMPath->primitive(2)->tdcTimeStamp(),
0613                                                                              latcomb[2],
0614                                                                              inMPath->primitive(3)->channelId(),
0615                                                                              inMPath->primitive(3)->tdcTimeStamp(),
0616                                                                              latcomb[3],
0617                                                                              -1}));
0618           }
0619         }
0620       }
0621     }
0622   } else {
0623     // implementation for advanced (2SL) grouping, no filter required..
0624     if (debug_)
0625       LogDebug("DTTrigPhase2Prod") << "Fitting 2SL at once ";
0626     for (auto& ch_muonpaths : muonpaths) {
0627       mpathanalyzer_->run(iEvent, iEventSetup, ch_muonpaths.second, outmpaths[ch_muonpaths.first]);
0628     }
0629   }
0630 
0631   skip_processing_ = skip_processing_ || output_slfitter_;
0632 
0633   if (dump_) {
0634     for (auto& ch_outmpaths : outmpaths) {
0635       for (unsigned int i = 0; i < ch_outmpaths.second.size(); i++) {
0636         LogInfo("DTTrigPhase2Prod") << iEvent.id().event() << " mp " << i << ": "
0637                                     << ch_outmpaths.second.at(i)->bxTimeValue() << " "
0638                                     << ch_outmpaths.second.at(i)->horizPos() << " "
0639                                     << ch_outmpaths.second.at(i)->tanPhi() << " " << ch_outmpaths.second.at(i)->phi()
0640                                     << " " << ch_outmpaths.second.at(i)->phiB() << " "
0641                                     << ch_outmpaths.second.at(i)->quality() << " "
0642                                     << ch_outmpaths.second.at(i)->chiSquare();
0643       }
0644     }
0645     for (auto& ch_metaPrimitives : metaPrimitives) {
0646       for (unsigned int i = 0; i < ch_metaPrimitives.second.size(); i++) {
0647         stringstream ss;
0648         ss << iEvent.id().event() << " mp " << i << ": ";
0649         printmP(ss.str(), ch_metaPrimitives.second.at(i));
0650       }
0651     }
0652   }
0653 
0654   muonpaths.clear();
0655   filteredmuonpaths.clear();
0656 
0657   /////////////////////////////////////
0658   //// CONFIRMATION:
0659   /////////////////////////////////////
0660 
0661   std::map<int, std::vector<metaPrimitive>> confirmedMetaPrimitives;
0662   for (auto& ch_metaPrimitives : metaPrimitives) {
0663     if (!skip_processing_ && allow_confirmation_)
0664       mpathconfirmator_->run(
0665           iEvent, iEventSetup, ch_metaPrimitives.second, dtdigis, confirmedMetaPrimitives[ch_metaPrimitives.first]);
0666     else
0667       for (auto& mp : ch_metaPrimitives.second) {
0668         confirmedMetaPrimitives[ch_metaPrimitives.first].push_back(mp);
0669       }
0670   }
0671 
0672   metaPrimitives.clear();
0673   skip_processing_ = skip_processing_ || output_confirmed_;
0674 
0675   /////////////////////////////////////
0676   //  FILTER SECTIONS:
0677   ////////////////////////////////////
0678 
0679   if (debug_)
0680     LogDebug("DTTrigPhase2Prod") << "declaring new vector for filtered" << std::endl;
0681 
0682   std::map<int, std::vector<metaPrimitive>> filteredMetaPrimitives;
0683   if (algo_ == Standard)
0684     for (auto& ch_confirmedMetaPrimitives : confirmedMetaPrimitives) {
0685       if (!skip_processing_)
0686         mpathqualityenhancer_->run(iEvent,
0687                                    iEventSetup,
0688                                    ch_confirmedMetaPrimitives.second,
0689                                    filteredMetaPrimitives[ch_confirmedMetaPrimitives.first]);
0690       else
0691         for (auto& mp : ch_confirmedMetaPrimitives.second) {
0692           filteredMetaPrimitives[ch_confirmedMetaPrimitives.first].push_back(mp);
0693         }
0694     }
0695   if (dump_) {
0696     for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0697       for (unsigned int i = 0; i < ch_filteredMetaPrimitives.second.size(); i++) {
0698         stringstream ss;
0699         ss << iEvent.id().event() << " filtered mp " << i << ": ";
0700         printmP(ss.str(), ch_filteredMetaPrimitives.second.at(i));
0701       }
0702     }
0703   }
0704 
0705   skip_processing_ = skip_processing_ || output_slfilter_;
0706   confirmedMetaPrimitives.clear();
0707 
0708   if (debug_)
0709     for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0710       LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0711                                    << ch_filteredMetaPrimitives.second.size() << " filteredMetaPrimitives (superlayer)"
0712                                    << std::endl;
0713     }
0714   if (debug_)
0715     LogDebug("DTTrigPhase2Prod") << "filteredMetaPrimitives: starting correlations" << std::endl;
0716 
0717   /////////////////////////////////////
0718   //// CORRELATION:
0719   /////////////////////////////////////
0720 
0721   std::map<int, std::vector<metaPrimitive>> correlatedMetaPrimitives;
0722   if (algo_ == Standard) {
0723     for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0724       if (!skip_processing_)
0725         mpathassociator_->run(iEvent,
0726                               iEventSetup,
0727                               ch_filteredMetaPrimitives.second,
0728                               correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
0729       else
0730         for (auto& mp : ch_filteredMetaPrimitives.second) {
0731           correlatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0732         }
0733     }
0734   } else {
0735     for (auto& ch_outmpaths : outmpaths) {
0736       for (const auto& muonpath : ch_outmpaths.second) {
0737         correlatedMetaPrimitives[ch_outmpaths.first].emplace_back(muonpath->rawId(),
0738                                                                   (double)muonpath->bxTimeValue(),
0739                                                                   muonpath->horizPos(),
0740                                                                   muonpath->tanPhi(),
0741                                                                   muonpath->phi(),
0742                                                                   muonpath->phiB(),
0743                                                                   muonpath->phi_cmssw(),
0744                                                                   muonpath->phiB_cmssw(),
0745                                                                   muonpath->chiSquare(),
0746                                                                   (int)muonpath->quality(),
0747                                                                   muonpath->primitive(0)->channelId(),
0748                                                                   muonpath->primitive(0)->tdcTimeStamp(),
0749                                                                   muonpath->primitive(0)->laterality(),
0750                                                                   muonpath->primitive(1)->channelId(),
0751                                                                   muonpath->primitive(1)->tdcTimeStamp(),
0752                                                                   muonpath->primitive(1)->laterality(),
0753                                                                   muonpath->primitive(2)->channelId(),
0754                                                                   muonpath->primitive(2)->tdcTimeStamp(),
0755                                                                   muonpath->primitive(2)->laterality(),
0756                                                                   muonpath->primitive(3)->channelId(),
0757                                                                   muonpath->primitive(3)->tdcTimeStamp(),
0758                                                                   muonpath->primitive(3)->laterality(),
0759                                                                   muonpath->primitive(4)->channelId(),
0760                                                                   muonpath->primitive(4)->tdcTimeStamp(),
0761                                                                   muonpath->primitive(4)->laterality(),
0762                                                                   muonpath->primitive(5)->channelId(),
0763                                                                   muonpath->primitive(5)->tdcTimeStamp(),
0764                                                                   muonpath->primitive(5)->laterality(),
0765                                                                   muonpath->primitive(6)->channelId(),
0766                                                                   muonpath->primitive(6)->tdcTimeStamp(),
0767                                                                   muonpath->primitive(6)->laterality(),
0768                                                                   muonpath->primitive(7)->channelId(),
0769                                                                   muonpath->primitive(7)->tdcTimeStamp(),
0770                                                                   muonpath->primitive(7)->laterality());
0771       }
0772     }
0773   }
0774 
0775   skip_processing_ = skip_processing_ || output_matcher_;
0776 
0777   if (debug_)
0778     for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0779       LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0780                                    << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
0781     }
0782   if (dump_) {
0783     for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0784       LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0785                                    << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
0786     }
0787     for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0788       for (unsigned int i = 0; i < ch_correlatedMetaPrimitives.second.size(); i++) {
0789         stringstream ss;
0790         ss << iEvent.id().event() << " correlated mp " << i << ": ";
0791         printmPC(ss.str(), ch_correlatedMetaPrimitives.second.at(i));
0792       }
0793     }
0794   }
0795 
0796   // Correlated Filtering
0797   std::map<int, std::vector<metaPrimitive>> filtCorrelatedMetaPrimitives;
0798   if (algo_ == Standard) {
0799     for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0800       if (!skip_processing_)
0801         mpathcorfilter_->run(iEvent,
0802                              iEventSetup,
0803                              ch_filteredMetaPrimitives.second,
0804                              correlatedMetaPrimitives[ch_filteredMetaPrimitives.first],
0805                              filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
0806       else {
0807         for (auto& mp : ch_filteredMetaPrimitives.second) {
0808           filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0809         }
0810         if (output_matcher_)
0811           for (auto& mp : correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]) {
0812             filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0813           }
0814       }
0815     }
0816   }
0817 
0818   correlatedMetaPrimitives.clear();
0819   filteredMetaPrimitives.clear();
0820 
0821   double shift_back = 0;
0822   if (scenario_ == MC)  //scope for MC
0823     shift_back = 400;
0824   else if (scenario_ == DATA)  //scope for data
0825     shift_back = 0;
0826   else if (scenario_ == SLICE_TEST)  //scope for slice test
0827     shift_back = 400;
0828 
0829   // RPC integration
0830   if (useRPC_) {
0831     rpc_integrator_->initialise(iEventSetup, shift_back);
0832     rpc_integrator_->prepareMetaPrimitives(rpcRecHits);
0833     for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
0834       rpc_integrator_->matchWithDTAndUseRPCTime(ch_correlatedMetaPrimitives.second);  // Probably this is a FIXME
0835     }
0836     rpc_integrator_->makeRPCOnlySegments();
0837     rpc_integrator_->storeRPCSingleHits();
0838     rpc_integrator_->removeRPCHitsUsed();
0839   }
0840 
0841   /// STORING RESULTs
0842   vector<L1Phase2MuDTPhDigi> outP2Ph;
0843   vector<L1Phase2MuDTExtPhDigi> outExtP2Ph;
0844   vector<L1Phase2MuDTThDigi> outP2Th;
0845   vector<L1Phase2MuDTExtThDigi> outExtP2Th;
0846 
0847   // Assigning index value
0848   if (!skip_processing_)
0849     for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
0850       assignIndex(ch_correlatedMetaPrimitives.second);
0851     }
0852 
0853   for (auto& ch_correlatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
0854     for (const auto& metaPrimitiveIt : ch_correlatedMetaPrimitives.second) {
0855       DTChamberId chId(metaPrimitiveIt.rawId);
0856       DTSuperLayerId slId(metaPrimitiveIt.rawId);
0857       if (debug_)
0858         LogDebug("DTTrigPhase2Prod") << "looping in final vector: SuperLayerId" << chId << " x=" << metaPrimitiveIt.x
0859                                      << " quality=" << metaPrimitiveIt.quality
0860                                      << " BX=" << round(metaPrimitiveIt.t0 / 25.) << " index=" << metaPrimitiveIt.index;
0861 
0862       int sectorTP = chId.sector();
0863       //sectors 13 and 14 exist only for the outermost stations for sectors 4 and 10 respectively
0864       //due to the larger MB4 that are divided into two.
0865       if (sectorTP == 13)
0866         sectorTP = 4;
0867       if (sectorTP == 14)
0868         sectorTP = 10;
0869       sectorTP = sectorTP - 1;
0870       int sl = 0;
0871       if (metaPrimitiveIt.quality < LOWLOWQ || metaPrimitiveIt.quality == CHIGHQ) {
0872         if (inner(metaPrimitiveIt))
0873           sl = 1;
0874         else
0875           sl = 3;
0876       }
0877 
0878       float tp_t0 =
0879           (metaPrimitiveIt.t0 - shift_back * LHC_CLK_FREQ) * ((float)TIME_TO_TDC_COUNTS / (float)LHC_CLK_FREQ);
0880 
0881       if (debug_)
0882         LogDebug("DTTrigPhase2Prod") << "pushing back phase-2 dataformat carlo-federica dataformat";
0883 
0884       if (slId.superLayer() != 2) {
0885         if (df_extended_ == 1 || df_extended_ == 2) {
0886           int pathWireId[8] = {metaPrimitiveIt.wi1,
0887                                metaPrimitiveIt.wi2,
0888                                metaPrimitiveIt.wi3,
0889                                metaPrimitiveIt.wi4,
0890                                metaPrimitiveIt.wi5,
0891                                metaPrimitiveIt.wi6,
0892                                metaPrimitiveIt.wi7,
0893                                metaPrimitiveIt.wi8};
0894 
0895           int pathTDC[8] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
0896                             max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
0897                             max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
0898                             max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1),
0899                             max((int)round(metaPrimitiveIt.tdc5 - shift_back * LHC_CLK_FREQ), -1),
0900                             max((int)round(metaPrimitiveIt.tdc6 - shift_back * LHC_CLK_FREQ), -1),
0901                             max((int)round(metaPrimitiveIt.tdc7 - shift_back * LHC_CLK_FREQ), -1),
0902                             max((int)round(metaPrimitiveIt.tdc8 - shift_back * LHC_CLK_FREQ), -1)};
0903 
0904           int pathLat[8] = {metaPrimitiveIt.lat1,
0905                             metaPrimitiveIt.lat2,
0906                             metaPrimitiveIt.lat3,
0907                             metaPrimitiveIt.lat4,
0908                             metaPrimitiveIt.lat5,
0909                             metaPrimitiveIt.lat6,
0910                             metaPrimitiveIt.lat7,
0911                             metaPrimitiveIt.lat8};
0912 
0913           // phiTP (extended DF)
0914           outExtP2Ph.emplace_back(
0915               L1Phase2MuDTExtPhDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0916                                     chId.wheel(),                                           // uwh   (m_wheel)
0917                                     sectorTP,                                               // usc   (m_sector)
0918                                     chId.station(),                                         // ust   (m_station)
0919                                     sl,                                                     // ust   (m_station)
0920                                     (int)round(metaPrimitiveIt.phi * PHIRES_CONV),          // uphi  (m_phiAngle)
0921                                     (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV),        // uphib (m_phiBending)
0922                                     metaPrimitiveIt.quality,                                // uqua  (m_qualityCode)
0923                                     metaPrimitiveIt.index,                                  // uind  (m_segmentIndex)
0924                                     tp_t0,                                                  // ut0   (m_t0Segment)
0925                                     (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),        // uchi2 (m_chi2Segment)
0926                                     (int)round(metaPrimitiveIt.x * 1000),                   // ux    (m_xLocal)
0927                                     (int)round(metaPrimitiveIt.tanPhi * 1000),              // utan  (m_tanPsi)
0928                                     (int)round(metaPrimitiveIt.phi_cmssw * PHIRES_CONV),    // uphi  (m_phiAngleCMSSW)
0929                                     (int)round(metaPrimitiveIt.phiB_cmssw * PHIBRES_CONV),  // uphib (m_phiBendingCMSSW)
0930                                     metaPrimitiveIt.rpcFlag,                                // urpc  (m_rpcFlag)
0931                                     pathWireId,
0932                                     pathTDC,
0933                                     pathLat));
0934         }
0935         if (df_extended_ == 0 || df_extended_ == 2) {
0936           // phiTP (standard DF)
0937           outP2Ph.push_back(L1Phase2MuDTPhDigi(
0938               (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0939               chId.wheel(),                                     // uwh (m_wheel)
0940               sectorTP,                                         // usc (m_sector)
0941               chId.station(),                                   // ust (m_station)
0942               sl,                                               // ust (m_station)
0943               (int)round(metaPrimitiveIt.phi * PHIRES_CONV),    // uphi (_phiAngle)
0944               (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV),  // uphib (m_phiBending)
0945               metaPrimitiveIt.quality,                          // uqua (m_qualityCode)
0946               metaPrimitiveIt.index,                            // uind (m_segmentIndex)
0947               tp_t0,                                            // ut0 (m_t0Segment)
0948               (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),  // uchi2 (m_chi2Segment)
0949               metaPrimitiveIt.rpcFlag                           // urpc (m_rpcFlag)
0950               ));
0951         }
0952       } else {
0953         if (df_extended_ == 1 || df_extended_ == 2) {
0954           int pathWireId[4] = {metaPrimitiveIt.wi1, metaPrimitiveIt.wi2, metaPrimitiveIt.wi3, metaPrimitiveIt.wi4};
0955 
0956           int pathTDC[4] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
0957                             max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
0958                             max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
0959                             max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1)};
0960 
0961           int pathLat[4] = {metaPrimitiveIt.lat1, metaPrimitiveIt.lat2, metaPrimitiveIt.lat3, metaPrimitiveIt.lat4};
0962 
0963           // thTP (extended DF)
0964           outExtP2Th.emplace_back(
0965               L1Phase2MuDTExtThDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0966                                     chId.wheel(),                                        // uwh   (m_wheel)
0967                                     sectorTP,                                            // usc   (m_sector)
0968                                     chId.station(),                                      // ust   (m_station)
0969                                     (int)round(metaPrimitiveIt.phi * ZRES_CONV),         // uz    (m_zGlobal)
0970                                     (int)round(metaPrimitiveIt.phiB * KRES_CONV),        // uk    (m_kSlope)
0971                                     metaPrimitiveIt.quality,                             // uqua  (m_qualityCode)
0972                                     metaPrimitiveIt.index,                               // uind  (m_segmentIndex)
0973                                     tp_t0,                                               // ut0   (m_t0Segment)
0974                                     (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),     // uchi2 (m_chi2Segment)
0975                                     (int)round(metaPrimitiveIt.x * 1000),                // ux    (m_yLocal)
0976                                     (int)round(metaPrimitiveIt.phi_cmssw * ZRES_CONV),   // uphi  (m_zCMSSW)
0977                                     (int)round(metaPrimitiveIt.phiB_cmssw * KRES_CONV),  // uphib (m_kCMSSW)
0978                                     metaPrimitiveIt.rpcFlag,                             // urpc  (m_rpcFlag)
0979                                     pathWireId,
0980                                     pathTDC,
0981                                     pathLat));
0982         }
0983         if (df_extended_ == 0 || df_extended_ == 2) {
0984           // thTP (standard DF)
0985           outP2Th.push_back(L1Phase2MuDTThDigi(
0986               (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0987               chId.wheel(),                                     // uwh (m_wheel)
0988               sectorTP,                                         // usc (m_sector)
0989               chId.station(),                                   // ust (m_station)
0990               (int)round(metaPrimitiveIt.phi * ZRES_CONV),      // uz (m_zGlobal)
0991               (int)round(metaPrimitiveIt.phiB * KRES_CONV),     // uk (m_kSlope)
0992               metaPrimitiveIt.quality,                          // uqua (m_qualityCode)
0993               metaPrimitiveIt.index,                            // uind (m_segmentIndex)
0994               tp_t0,                                            // ut0 (m_t0Segment)
0995               (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),  // uchi2 (m_chi2Segment)
0996               metaPrimitiveIt.rpcFlag                           // urpc (m_rpcFlag)
0997               ));
0998         }
0999       }
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 
1056 bool DTTrigPhase2Prod::outer(const metaPrimitive& mp) const {
1057   int counter = (mp.wi5 != -1) + (mp.wi6 != -1) + (mp.wi7 != -1) + (mp.wi8 != -1);
1058   return (counter > 2);
1059 }
1060 
1061 bool DTTrigPhase2Prod::inner(const metaPrimitive& mp) const {
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 
1078 void DTTrigPhase2Prod::printmP(const metaPrimitive& mP) const {
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 
1105 void DTTrigPhase2Prod::printmPC(const metaPrimitive& mP) const {
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 
1122 int DTTrigPhase2Prod::rango(const metaPrimitive& mp) const {
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 
1178 int DTTrigPhase2Prod::assignQualityOrder(const metaPrimitive& mP) const {
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 
1234 void DTTrigPhase2Prod::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
1235   // dtTriggerPhase2PrimitiveDigis
1236   edm::ParameterSetDescription desc;
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   {
1292     edm::ParameterSetDescription psd0;
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   {
1310     edm::ParameterSetDescription psd0;
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 
1329 DEFINE_FWK_MODULE(DTTrigPhase2Prod);