File indexing completed on 2024-09-07 04:36: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/MPCoincidenceFilter.h"
0043 #include "L1Trigger/DTTriggerPhase2/interface/MPQualityEnhancerFilter.h"
0044 #include "L1Trigger/DTTriggerPhase2/interface/MPRedundantFilter.h"
0045 #include "L1Trigger/DTTriggerPhase2/interface/MPCleanHitsFilter.h"
0046 #include "L1Trigger/DTTriggerPhase2/interface/MPQualityEnhancerFilterBayes.h"
0047 #include "L1Trigger/DTTriggerPhase2/interface/GlobalCoordsObtainer.h"
0048
0049 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0050 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
0051 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0052 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0053 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
0054 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTPhContainer.h"
0055 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTPhDigi.h"
0056 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtPhContainer.h"
0057 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtPhDigi.h"
0058 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTThContainer.h"
0059 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTThDigi.h"
0060 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtThContainer.h"
0061 #include "DataFormats/L1DTTrackFinder/interface/L1Phase2MuDTExtThDigi.h"
0062
0063
0064 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
0065
0066
0067 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
0068 #include <DataFormats/MuonDetId/interface/RPCDetId.h>
0069 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0070 #include "L1Trigger/DTTriggerPhase2/interface/RPCIntegrator.h"
0071
0072 #include <fstream>
0073 #include <iostream>
0074 #include <queue>
0075 #include <cmath>
0076
0077 using namespace edm;
0078 using namespace std;
0079 using namespace cmsdt;
0080
0081 class DTTrigPhase2Prod : public edm::stream::EDProducer<> {
0082 typedef std::map<DTChamberId, DTDigiCollection, std::less<DTChamberId>> DTDigiMap;
0083 typedef DTDigiMap::iterator DTDigiMap_iterator;
0084 typedef DTDigiMap::const_iterator DTDigiMap_const_iterator;
0085
0086 public:
0087
0088 DTTrigPhase2Prod(const edm::ParameterSet& pset);
0089
0090
0091 ~DTTrigPhase2Prod() override;
0092
0093
0094 void beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
0095
0096
0097 void produce(edm::Event& iEvent, const edm::EventSetup& iEventSetup) override;
0098
0099
0100 void endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
0101
0102
0103 int rango(const metaPrimitive& mp) const;
0104 bool outer(const metaPrimitive& mp) const;
0105 bool inner(const metaPrimitive& mp) const;
0106 void printmP(const std::string& ss, const metaPrimitive& mP) const;
0107 void printmP(const metaPrimitive& mP) const;
0108 void printmPC(const std::string& ss, const metaPrimitive& mP) const;
0109 void printmPC(const metaPrimitive& mP) const;
0110 bool hasPosRF(int wh, int sec) const;
0111
0112
0113 MP_QUALITY getMinimumQuality(void);
0114
0115
0116 void setChiSquareThreshold(float ch2Thr);
0117 void setMinimumQuality(MP_QUALITY q);
0118
0119 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0120
0121
0122 const DTGeometry* dtGeo_;
0123 edm::ESGetToken<DTGeometry, MuonGeometryRecord> dtGeomH;
0124 std::vector<std::pair<int, MuonPath>> primitives_;
0125
0126 private:
0127
0128 bool my_CCBValid_;
0129
0130
0131 int my_BXoffset_;
0132
0133
0134 bool debug_;
0135 bool dump_;
0136 double dT0_correlate_TP_;
0137 int scenario_;
0138 int df_extended_;
0139 int co_option_;
0140 int co_quality_;
0141 int max_index_;
0142
0143 bool output_mixer_;
0144 bool output_latpredictor_;
0145 bool output_slfitter_;
0146 bool output_slfilter_;
0147 bool output_confirmed_;
0148 bool output_matcher_;
0149 bool skip_processing_;
0150 bool allow_confirmation_;
0151
0152
0153 edm::EDGetTokenT<DTDigiCollection> dtDigisToken_;
0154 edm::EDGetTokenT<RPCRecHitCollection> rpcRecHitsLabel_;
0155
0156
0157 int algo_;
0158 std::unique_ptr<MotherGrouping> grouping_obj_;
0159 std::unique_ptr<MuonPathAnalyzer> mpathanalyzer_;
0160 std::unique_ptr<LateralityProvider> latprovider_;
0161 std::unique_ptr<MPFilter> mpathqualityenhancer_;
0162 std::unique_ptr<MPFilter> mpathqualityenhancerbayes_;
0163 std::unique_ptr<MPFilter> mpathredundantfilter_;
0164 std::unique_ptr<MPFilter> mpathhitsfilter_;
0165 std::unique_ptr<MuonPathAnalyzer> mpathassociator_;
0166 std::unique_ptr<MuonPathConfirmator> mpathconfirmator_;
0167 std::unique_ptr<MPFilter> mpathcorfilter_;
0168 std::unique_ptr<MPFilter> mpathcoifilter_;
0169 std::shared_ptr<GlobalCoordsObtainer> globalcoordsobtainer_;
0170
0171
0172 bool activateBuffer_;
0173 int superCellhalfspacewidth_;
0174 float superCelltimewidth_;
0175 std::vector<DTDigiCollection*> distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ);
0176 void processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
0177 std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec);
0178
0179
0180 std::unique_ptr<RPCIntegrator> rpc_integrator_;
0181 bool useRPC_;
0182
0183 void assignIndex(std::vector<metaPrimitive>& inMPaths);
0184 void assignIndexPerBX(std::vector<metaPrimitive>& inMPaths);
0185 int assignQualityOrder(const metaPrimitive& mP) const;
0186
0187 const std::unordered_map<int, int> qmap_;
0188 };
0189
0190 namespace {
0191 struct {
0192 bool operator()(std::pair<DTLayerId, DTDigi> a, std::pair<DTLayerId, DTDigi> b) const {
0193 return (a.second.time() < b.second.time());
0194 }
0195 } const DigiTimeOrdering;
0196 }
0197
0198 DTTrigPhase2Prod::DTTrigPhase2Prod(const ParameterSet& pset)
0199 : qmap_({{8, 8}, {7, 7}, {6, 6}, {4, 4}, {3, 3}, {2, 2}, {1, 1}}) {
0200 produces<L1Phase2MuDTPhContainer>();
0201 produces<L1Phase2MuDTThContainer>();
0202 produces<L1Phase2MuDTExtPhContainer>();
0203 produces<L1Phase2MuDTExtThContainer>();
0204
0205 debug_ = pset.getUntrackedParameter<bool>("debug");
0206 dump_ = pset.getUntrackedParameter<bool>("dump");
0207
0208 scenario_ = pset.getParameter<int>("scenario");
0209
0210 df_extended_ = pset.getParameter<int>("df_extended");
0211 co_option_ = pset.getParameter<int>("co_option");
0212 co_quality_ = pset.getParameter<int>("co_quality");
0213 max_index_ = pset.getParameter<int>("max_primitives") - 1;
0214
0215 dtDigisToken_ = consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiTag"));
0216
0217 rpcRecHitsLabel_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("rpcRecHits"));
0218 useRPC_ = pset.getParameter<bool>("useRPC");
0219
0220
0221 algo_ = pset.getParameter<int>("algo");
0222
0223
0224
0225 output_mixer_ = pset.getParameter<bool>("output_mixer");
0226 output_latpredictor_ = pset.getParameter<bool>("output_latpredictor");
0227 output_slfitter_ = pset.getParameter<bool>("output_slfitter");
0228 output_slfilter_ = pset.getParameter<bool>("output_slfilter");
0229 output_confirmed_ = pset.getParameter<bool>("output_confirmed");
0230 output_matcher_ = pset.getParameter<bool>("output_matcher");
0231 allow_confirmation_ = pset.getParameter<bool>("allow_confirmation");
0232
0233 edm::ConsumesCollector consumesColl(consumesCollector());
0234 globalcoordsobtainer_ = std::make_shared<GlobalCoordsObtainer>(pset);
0235 globalcoordsobtainer_->generate_luts();
0236
0237 if (algo_ == PseudoBayes) {
0238 grouping_obj_ =
0239 std::make_unique<PseudoBayesGrouping>(pset.getParameter<edm::ParameterSet>("PseudoBayesPattern"), consumesColl);
0240 } else if (algo_ == HoughTrans) {
0241 grouping_obj_ =
0242 std::make_unique<HoughGrouping>(pset.getParameter<edm::ParameterSet>("HoughGrouping"), consumesColl);
0243 } else {
0244 grouping_obj_ = std::make_unique<TrapezoidalGrouping>(pset, consumesColl);
0245 }
0246
0247 if (algo_ == Standard) {
0248 if (debug_)
0249 LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: JM analyzer";
0250 mpathanalyzer_ = std::make_unique<MuonPathSLFitter>(pset, consumesColl, globalcoordsobtainer_);
0251 latprovider_ = std::make_unique<LateralityCoarsedProvider>(pset, consumesColl);
0252 } else {
0253 if (debug_)
0254 LogDebug("DTTrigPhase2Prod") << "DTp2:constructor: Full chamber analyzer";
0255 mpathanalyzer_ = std::make_unique<MuonPathAnalyzerInChamber>(pset, consumesColl, globalcoordsobtainer_);
0256 }
0257
0258
0259 activateBuffer_ = pset.getParameter<bool>("activateBuffer");
0260 superCellhalfspacewidth_ = pset.getParameter<int>("superCellspacewidth") / 2;
0261 superCelltimewidth_ = pset.getParameter<double>("superCelltimewidth");
0262
0263 mpathqualityenhancer_ = std::make_unique<MPSLFilter>(pset);
0264 mpathqualityenhancerbayes_ = std::make_unique<MPQualityEnhancerFilterBayes>(pset);
0265 mpathredundantfilter_ = std::make_unique<MPRedundantFilter>(pset);
0266 mpathhitsfilter_ = std::make_unique<MPCleanHitsFilter>(pset);
0267 mpathconfirmator_ = std::make_unique<MuonPathConfirmator>(pset, consumesColl);
0268 mpathassociator_ = std::make_unique<MuonPathCorFitter>(pset, consumesColl, globalcoordsobtainer_);
0269 mpathcorfilter_ = std::make_unique<MPCorFilter>(pset);
0270 mpathcoifilter_ = std::make_unique<MPCoincidenceFilter>(pset);
0271 rpc_integrator_ = std::make_unique<RPCIntegrator>(pset, consumesColl);
0272
0273 dtGeomH = esConsumes<DTGeometry, MuonGeometryRecord, edm::Transition::BeginRun>();
0274 }
0275
0276 DTTrigPhase2Prod::~DTTrigPhase2Prod() {
0277 if (debug_)
0278 LogDebug("DTTrigPhase2Prod") << "DTp2: calling destructor" << std::endl;
0279 }
0280
0281 void DTTrigPhase2Prod::beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
0282 if (debug_)
0283 LogDebug("DTTrigPhase2Prod") << "beginRun " << iRun.id().run();
0284 if (debug_)
0285 LogDebug("DTTrigPhase2Prod") << "beginRun: getting DT geometry";
0286
0287 grouping_obj_->initialise(iEventSetup);
0288 mpathanalyzer_->initialise(iEventSetup);
0289 mpathqualityenhancer_->initialise(iEventSetup);
0290 mpathredundantfilter_->initialise(iEventSetup);
0291 mpathqualityenhancerbayes_->initialise(iEventSetup);
0292 mpathhitsfilter_->initialise(iEventSetup);
0293 mpathassociator_->initialise(iEventSetup);
0294 mpathcorfilter_->initialise(iEventSetup);
0295 mpathcoifilter_->initialise(iEventSetup);
0296
0297 if (auto geom = iEventSetup.getHandle(dtGeomH)) {
0298 dtGeo_ = &(*geom);
0299 }
0300 }
0301
0302 void DTTrigPhase2Prod::produce(Event& iEvent, const EventSetup& iEventSetup) {
0303 if (debug_)
0304 LogDebug("DTTrigPhase2Prod") << "produce";
0305 edm::Handle<DTDigiCollection> dtdigis;
0306 iEvent.getByToken(dtDigisToken_, dtdigis);
0307
0308 if (debug_)
0309 LogDebug("DTTrigPhase2Prod") << "\t Getting the RPC RecHits" << std::endl;
0310 edm::Handle<RPCRecHitCollection> rpcRecHits;
0311 iEvent.getByToken(rpcRecHitsLabel_, rpcRecHits);
0312
0313
0314
0315
0316
0317 DTDigiMap digiMap;
0318 DTDigiCollection::DigiRangeIterator detUnitIt;
0319 for (const auto& detUnitIt : *dtdigis) {
0320 const DTLayerId& layId = detUnitIt.first;
0321 const DTChamberId chambId = layId.superlayerId().chamberId();
0322 const DTDigiCollection::Range& range = detUnitIt.second;
0323 digiMap[chambId].put(range, layId);
0324 }
0325
0326
0327 if (debug_ && activateBuffer_)
0328 LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber using a buffer and super cells.";
0329 else if (debug_)
0330 LogDebug("DTTrigPhase2Prod") << "produce - Getting and grouping digis per chamber.";
0331
0332 std::map<int, MuonPathPtrs> muonpaths;
0333 for (const auto& ich : dtGeo_->chambers()) {
0334
0335 const DTChamber* chamb = ich;
0336 DTChamberId chid = chamb->id();
0337 DTDigiMap_iterator dmit = digiMap.find(chid);
0338
0339 if (dmit == digiMap.end())
0340 continue;
0341
0342 if (activateBuffer_) {
0343
0344 std::vector<std::pair<DTLayerId, DTDigi>> tmpvec;
0345 tmpvec.clear();
0346
0347 for (const auto& dtLayerIdIt : (*dmit).second) {
0348 for (DTDigiCollection::const_iterator digiIt = (dtLayerIdIt.second).first;
0349 digiIt != (dtLayerIdIt.second).second;
0350 digiIt++) {
0351 tmpvec.emplace_back(dtLayerIdIt.first, *digiIt);
0352 }
0353 }
0354
0355
0356 if (tmpvec.empty())
0357 continue;
0358
0359
0360 std::sort(tmpvec.begin(), tmpvec.end(), DigiTimeOrdering);
0361 std::queue<std::pair<DTLayerId, DTDigi>> timequeue;
0362
0363 for (const auto& elem : tmpvec)
0364 timequeue.emplace(elem);
0365 tmpvec.clear();
0366
0367
0368 std::vector<DTDigiCollection*> superCells;
0369 superCells = distribDigis(timequeue);
0370
0371
0372
0373
0374 while (!superCells.empty()) {
0375 grouping_obj_->run(iEvent, iEventSetup, *(superCells.back()), muonpaths[chid.rawId()]);
0376 superCells.pop_back();
0377 }
0378 } else {
0379 grouping_obj_->run(iEvent, iEventSetup, (*dmit).second, muonpaths[chid.rawId()]);
0380 }
0381 }
0382 digiMap.clear();
0383
0384 if (dump_) {
0385 for (auto& ch_muonpaths : muonpaths) {
0386 for (unsigned int i = 0; i < ch_muonpaths.second.size(); i++) {
0387 stringstream ss;
0388 ss << iEvent.id().event() << " mpath " << i << ": ";
0389 for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0390 ss << ch_muonpaths.second.at(i)->primitive(lay)->channelId() << " ";
0391 for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0392 ss << ch_muonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
0393 for (int lay = 0; lay < ch_muonpaths.second.at(i)->nprimitives(); lay++)
0394 ss << ch_muonpaths.second.at(i)->primitive(lay)->laterality() << " ";
0395 LogInfo("DTTrigPhase2Prod") << ss.str();
0396 }
0397 }
0398 }
0399
0400 std::map<int, std::vector<lat_vector>> lateralities;
0401 if (!output_mixer_) {
0402 for (auto& ch_muonpaths : muonpaths) {
0403 if (algo_ == Standard) {
0404 latprovider_->run(iEvent, iEventSetup, ch_muonpaths.second, lateralities[ch_muonpaths.first]);
0405 }
0406 }
0407 }
0408
0409
0410 std::map<int, MuonPathPtrs> filteredmuonpaths;
0411 for (auto& ch_muonpaths : muonpaths) {
0412 if (algo_ == Standard) {
0413 mpathredundantfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
0414 } else {
0415 mpathhitsfilter_->run(iEvent, iEventSetup, ch_muonpaths.second, filteredmuonpaths[ch_muonpaths.first]);
0416 }
0417 }
0418
0419 if (dump_) {
0420 for (auto& ch_filteredmuonpaths : filteredmuonpaths) {
0421 for (unsigned int i = 0; i < ch_filteredmuonpaths.second.size(); i++) {
0422 stringstream ss;
0423 ss << iEvent.id().event() << " filt. mpath " << i << ": ";
0424 for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
0425 ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->channelId() << " ";
0426 for (int lay = 0; lay < ch_filteredmuonpaths.second.at(i)->nprimitives(); lay++)
0427 ss << ch_filteredmuonpaths.second.at(i)->primitive(lay)->tdcTimeStamp() << " ";
0428 LogInfo("DTTrigPhase2Prod") << ss.str();
0429 }
0430 }
0431 }
0432
0433 skip_processing_ = output_mixer_ || output_latpredictor_;
0434
0435
0436
0437
0438
0439 if (debug_) {
0440 for (auto& ch_muonpaths : muonpaths) {
0441 LogDebug("DTTrigPhase2Prod") << "MUON PATHS found: " << ch_muonpaths.second.size() << " ("
0442 << filteredmuonpaths[ch_muonpaths.first].size() << ") in event "
0443 << iEvent.id().event();
0444 }
0445 }
0446 if (debug_)
0447 LogDebug("DTTrigPhase2Prod") << "filling NmetaPrimtives" << std::endl;
0448 std::map<int, std::vector<metaPrimitive>> metaPrimitives;
0449 std::map<int, MuonPathPtrs> outmpaths;
0450 if (algo_ == Standard) {
0451 if (debug_)
0452 LogDebug("DTTrigPhase2Prod") << "Fitting 1SL ";
0453 for (auto& ch_muonpaths : muonpaths) {
0454 if (!output_mixer_ && !output_latpredictor_)
0455 mpathanalyzer_->run(iEvent,
0456 iEventSetup,
0457 ch_muonpaths.second,
0458 lateralities[ch_muonpaths.first],
0459 metaPrimitives[ch_muonpaths.first]);
0460 else if (output_mixer_) {
0461 for (auto& inMPath : ch_muonpaths.second) {
0462 auto sl = inMPath->primitive(0)->superLayerId();
0463 int selected_lay = 1;
0464 if (inMPath->primitive(0)->tdcTimeStamp() != -1)
0465 selected_lay = 0;
0466 int dumLayId = inMPath->primitive(selected_lay)->cameraId();
0467 auto dtDumlayerId = DTLayerId(dumLayId);
0468 DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
0469 if (sl == 0)
0470 metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0471 -1,
0472 -1,
0473 -1,
0474 -1,
0475 -1,
0476 -1,
0477 -1,
0478 -1,
0479 -1,
0480 inMPath->primitive(0)->channelId(),
0481 inMPath->primitive(0)->tdcTimeStamp(),
0482 -1,
0483 inMPath->primitive(1)->channelId(),
0484 inMPath->primitive(1)->tdcTimeStamp(),
0485 -1,
0486 inMPath->primitive(2)->channelId(),
0487 inMPath->primitive(2)->tdcTimeStamp(),
0488 -1,
0489 inMPath->primitive(3)->channelId(),
0490 inMPath->primitive(3)->tdcTimeStamp(),
0491 -1,
0492 -1,
0493 -1,
0494 -1,
0495 -1,
0496 -1,
0497 -1,
0498 -1,
0499 -1,
0500 -1,
0501 -1,
0502 -1,
0503 -1,
0504 -1}));
0505 else
0506 metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
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 -1,
0521 -1,
0522 -1,
0523 -1,
0524 -1,
0525 -1,
0526 -1,
0527 -1,
0528 inMPath->primitive(0)->channelId(),
0529 inMPath->primitive(0)->tdcTimeStamp(),
0530 -1,
0531 inMPath->primitive(1)->channelId(),
0532 inMPath->primitive(1)->tdcTimeStamp(),
0533 -1,
0534 inMPath->primitive(2)->channelId(),
0535 inMPath->primitive(2)->tdcTimeStamp(),
0536 -1,
0537 inMPath->primitive(3)->channelId(),
0538 inMPath->primitive(3)->tdcTimeStamp(),
0539 -1,
0540 -1}));
0541 }
0542 } else if (output_latpredictor_) {
0543 int imp = -1;
0544 for (auto& inMPath : ch_muonpaths.second) {
0545 imp++;
0546 auto sl = inMPath->primitive(0)->superLayerId();
0547 int selected_lay = 1;
0548 if (inMPath->primitive(0)->tdcTimeStamp() != -1)
0549 selected_lay = 0;
0550 int dumLayId = inMPath->primitive(selected_lay)->cameraId();
0551 auto dtDumlayerId = DTLayerId(dumLayId);
0552 DTSuperLayerId MuonPathSLId(dtDumlayerId.wheel(), dtDumlayerId.station(), dtDumlayerId.sector(), sl + 1);
0553 for (auto& latcomb : lateralities[ch_muonpaths.first][imp]) {
0554 if (sl == 0)
0555 metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
0556 -1,
0557 -1,
0558 -1,
0559 -1,
0560 -1,
0561 -1,
0562 -1,
0563 -1,
0564 -1,
0565 inMPath->primitive(0)->channelId(),
0566 inMPath->primitive(0)->tdcTimeStamp(),
0567 latcomb[0],
0568 inMPath->primitive(1)->channelId(),
0569 inMPath->primitive(1)->tdcTimeStamp(),
0570 latcomb[1],
0571 inMPath->primitive(2)->channelId(),
0572 inMPath->primitive(2)->tdcTimeStamp(),
0573 latcomb[2],
0574 inMPath->primitive(3)->channelId(),
0575 inMPath->primitive(3)->tdcTimeStamp(),
0576 latcomb[3],
0577 -1,
0578 -1,
0579 -1,
0580 -1,
0581 -1,
0582 -1,
0583 -1,
0584 -1,
0585 -1,
0586 -1,
0587 -1,
0588 -1,
0589 -1}));
0590 else
0591 metaPrimitives[ch_muonpaths.first].emplace_back(metaPrimitive({MuonPathSLId.rawId(),
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 -1,
0606 -1,
0607 -1,
0608 -1,
0609 -1,
0610 -1,
0611 -1,
0612 -1,
0613 inMPath->primitive(0)->channelId(),
0614 inMPath->primitive(0)->tdcTimeStamp(),
0615 latcomb[0],
0616 inMPath->primitive(1)->channelId(),
0617 inMPath->primitive(1)->tdcTimeStamp(),
0618 latcomb[1],
0619 inMPath->primitive(2)->channelId(),
0620 inMPath->primitive(2)->tdcTimeStamp(),
0621 latcomb[2],
0622 inMPath->primitive(3)->channelId(),
0623 inMPath->primitive(3)->tdcTimeStamp(),
0624 latcomb[3],
0625 -1}));
0626 }
0627 }
0628 }
0629 }
0630 } else {
0631
0632 if (debug_)
0633 LogDebug("DTTrigPhase2Prod") << "Fitting 2SL at once ";
0634 for (auto& ch_muonpaths : muonpaths) {
0635 mpathanalyzer_->run(iEvent, iEventSetup, ch_muonpaths.second, outmpaths[ch_muonpaths.first]);
0636 }
0637 }
0638
0639 skip_processing_ = skip_processing_ || output_slfitter_;
0640
0641 if (dump_) {
0642 for (auto& ch_outmpaths : outmpaths) {
0643 for (unsigned int i = 0; i < ch_outmpaths.second.size(); i++) {
0644 LogInfo("DTTrigPhase2Prod") << iEvent.id().event() << " mp " << i << ": "
0645 << ch_outmpaths.second.at(i)->bxTimeValue() << " "
0646 << ch_outmpaths.second.at(i)->horizPos() << " "
0647 << ch_outmpaths.second.at(i)->tanPhi() << " " << ch_outmpaths.second.at(i)->phi()
0648 << " " << ch_outmpaths.second.at(i)->phiB() << " "
0649 << ch_outmpaths.second.at(i)->quality() << " "
0650 << ch_outmpaths.second.at(i)->chiSquare();
0651 }
0652 }
0653 for (auto& ch_metaPrimitives : metaPrimitives) {
0654 for (unsigned int i = 0; i < ch_metaPrimitives.second.size(); i++) {
0655 stringstream ss;
0656 ss << iEvent.id().event() << " mp " << i << ": ";
0657 printmP(ss.str(), ch_metaPrimitives.second.at(i));
0658 }
0659 }
0660 }
0661
0662 muonpaths.clear();
0663 filteredmuonpaths.clear();
0664
0665
0666
0667
0668
0669 std::map<int, std::vector<metaPrimitive>> confirmedMetaPrimitives;
0670 for (auto& ch_metaPrimitives : metaPrimitives) {
0671 if (!skip_processing_ && allow_confirmation_)
0672 mpathconfirmator_->run(
0673 iEvent, iEventSetup, ch_metaPrimitives.second, dtdigis, confirmedMetaPrimitives[ch_metaPrimitives.first]);
0674 else
0675 for (auto& mp : ch_metaPrimitives.second) {
0676 confirmedMetaPrimitives[ch_metaPrimitives.first].push_back(mp);
0677 }
0678 }
0679
0680 metaPrimitives.clear();
0681 skip_processing_ = skip_processing_ || output_confirmed_;
0682
0683
0684
0685
0686
0687 if (debug_)
0688 LogDebug("DTTrigPhase2Prod") << "declaring new vector for filtered" << std::endl;
0689
0690 std::map<int, std::vector<metaPrimitive>> filteredMetaPrimitives;
0691 if (algo_ == Standard)
0692 for (auto& ch_confirmedMetaPrimitives : confirmedMetaPrimitives) {
0693 if (!skip_processing_)
0694 mpathqualityenhancer_->run(iEvent,
0695 iEventSetup,
0696 ch_confirmedMetaPrimitives.second,
0697 filteredMetaPrimitives[ch_confirmedMetaPrimitives.first]);
0698 else
0699 for (auto& mp : ch_confirmedMetaPrimitives.second) {
0700 filteredMetaPrimitives[ch_confirmedMetaPrimitives.first].push_back(mp);
0701 }
0702 }
0703 if (dump_) {
0704 for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0705 for (unsigned int i = 0; i < ch_filteredMetaPrimitives.second.size(); i++) {
0706 stringstream ss;
0707 ss << iEvent.id().event() << " filtered mp " << i << ": ";
0708 printmP(ss.str(), ch_filteredMetaPrimitives.second.at(i));
0709 }
0710 }
0711 }
0712
0713 skip_processing_ = skip_processing_ || output_slfilter_;
0714 confirmedMetaPrimitives.clear();
0715
0716 if (debug_)
0717 for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0718 LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0719 << ch_filteredMetaPrimitives.second.size() << " filteredMetaPrimitives (superlayer)"
0720 << std::endl;
0721 }
0722 if (debug_)
0723 LogDebug("DTTrigPhase2Prod") << "filteredMetaPrimitives: starting correlations" << std::endl;
0724
0725
0726
0727
0728
0729 std::map<int, std::vector<metaPrimitive>> correlatedMetaPrimitives;
0730 if (algo_ == Standard) {
0731 for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0732 if (!skip_processing_)
0733 mpathassociator_->run(iEvent,
0734 iEventSetup,
0735 ch_filteredMetaPrimitives.second,
0736 correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
0737 else
0738 for (auto& mp : ch_filteredMetaPrimitives.second) {
0739 correlatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0740 }
0741 }
0742 } else {
0743 for (auto& ch_outmpaths : outmpaths) {
0744 for (const auto& muonpath : ch_outmpaths.second) {
0745 correlatedMetaPrimitives[ch_outmpaths.first].emplace_back(muonpath->rawId(),
0746 (double)muonpath->bxTimeValue(),
0747 muonpath->horizPos(),
0748 muonpath->tanPhi(),
0749 muonpath->phi(),
0750 muonpath->phiB(),
0751 muonpath->phi_cmssw(),
0752 muonpath->phiB_cmssw(),
0753 muonpath->chiSquare(),
0754 (int)muonpath->quality(),
0755 muonpath->primitive(0)->channelId(),
0756 muonpath->primitive(0)->tdcTimeStamp(),
0757 muonpath->primitive(0)->laterality(),
0758 muonpath->primitive(1)->channelId(),
0759 muonpath->primitive(1)->tdcTimeStamp(),
0760 muonpath->primitive(1)->laterality(),
0761 muonpath->primitive(2)->channelId(),
0762 muonpath->primitive(2)->tdcTimeStamp(),
0763 muonpath->primitive(2)->laterality(),
0764 muonpath->primitive(3)->channelId(),
0765 muonpath->primitive(3)->tdcTimeStamp(),
0766 muonpath->primitive(3)->laterality(),
0767 muonpath->primitive(4)->channelId(),
0768 muonpath->primitive(4)->tdcTimeStamp(),
0769 muonpath->primitive(4)->laterality(),
0770 muonpath->primitive(5)->channelId(),
0771 muonpath->primitive(5)->tdcTimeStamp(),
0772 muonpath->primitive(5)->laterality(),
0773 muonpath->primitive(6)->channelId(),
0774 muonpath->primitive(6)->tdcTimeStamp(),
0775 muonpath->primitive(6)->laterality(),
0776 muonpath->primitive(7)->channelId(),
0777 muonpath->primitive(7)->tdcTimeStamp(),
0778 muonpath->primitive(7)->laterality());
0779 }
0780 }
0781 }
0782
0783 skip_processing_ = skip_processing_ || output_matcher_;
0784
0785 if (debug_)
0786 for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0787 LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0788 << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
0789 }
0790 if (dump_) {
0791 for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0792 LogDebug("DTTrigPhase2Prod") << "DTp2 in event:" << iEvent.id().event() << " we found "
0793 << ch_correlatedMetaPrimitives.second.size() << " correlatedMetPrimitives (chamber)";
0794 }
0795 for (auto& ch_correlatedMetaPrimitives : correlatedMetaPrimitives) {
0796 for (unsigned int i = 0; i < ch_correlatedMetaPrimitives.second.size(); i++) {
0797 stringstream ss;
0798 ss << iEvent.id().event() << " correlated mp " << i << ": ";
0799 printmPC(ss.str(), ch_correlatedMetaPrimitives.second.at(i));
0800 }
0801 }
0802 }
0803
0804
0805 std::map<int, std::vector<metaPrimitive>> filtCorrelatedMetaPrimitives;
0806 if (algo_ == Standard) {
0807 for (auto& ch_filteredMetaPrimitives : filteredMetaPrimitives) {
0808 if (!skip_processing_)
0809 mpathcorfilter_->run(iEvent,
0810 iEventSetup,
0811 ch_filteredMetaPrimitives.second,
0812 correlatedMetaPrimitives[ch_filteredMetaPrimitives.first],
0813 filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first]);
0814 else {
0815 for (auto& mp : ch_filteredMetaPrimitives.second) {
0816 filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0817 }
0818 if (output_matcher_)
0819 for (auto& mp : correlatedMetaPrimitives[ch_filteredMetaPrimitives.first]) {
0820 filtCorrelatedMetaPrimitives[ch_filteredMetaPrimitives.first].push_back(mp);
0821 }
0822 }
0823 }
0824 }
0825
0826 correlatedMetaPrimitives.clear();
0827 filteredMetaPrimitives.clear();
0828
0829
0830
0831 std::vector<metaPrimitive> allMetaPrimitives;
0832 for (auto& ch_filtcorrelatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
0833 for (const auto& metaPrimitiveIt : ch_filtcorrelatedMetaPrimitives.second) {
0834 allMetaPrimitives.push_back(metaPrimitiveIt);
0835 }
0836 }
0837
0838 std::map<int, std::vector<metaPrimitive>> coMetaPrimitives;
0839 if (algo_ == Standard) {
0840 for (auto& ch_filtcorrelatedMetaPrimitives : filtCorrelatedMetaPrimitives) {
0841 if (!skip_processing_)
0842 mpathcoifilter_->run(iEvent,
0843 iEventSetup,
0844 allMetaPrimitives,
0845 filtCorrelatedMetaPrimitives[ch_filtcorrelatedMetaPrimitives.first],
0846 coMetaPrimitives[ch_filtcorrelatedMetaPrimitives.first]);
0847 else {
0848 for (auto& mp : ch_filtcorrelatedMetaPrimitives.second) {
0849 coMetaPrimitives[ch_filtcorrelatedMetaPrimitives.first].push_back(mp);
0850 }
0851 }
0852 }
0853 }
0854
0855 allMetaPrimitives.clear();
0856
0857
0858
0859 double shift_back = 0;
0860 if (scenario_ == MC)
0861 shift_back = 400;
0862 else if (scenario_ == DATA)
0863 shift_back = 0;
0864 else if (scenario_ == SLICE_TEST)
0865 shift_back = 400;
0866
0867
0868 if (useRPC_) {
0869 rpc_integrator_->initialise(iEventSetup, shift_back);
0870 rpc_integrator_->prepareMetaPrimitives(rpcRecHits);
0871 for (auto& ch_correlatedMetaPrimitives : coMetaPrimitives) {
0872 rpc_integrator_->matchWithDTAndUseRPCTime(ch_correlatedMetaPrimitives.second);
0873 }
0874 rpc_integrator_->makeRPCOnlySegments();
0875 rpc_integrator_->storeRPCSingleHits();
0876 rpc_integrator_->removeRPCHitsUsed();
0877 }
0878
0879
0880 vector<L1Phase2MuDTPhDigi> outP2Ph;
0881 vector<L1Phase2MuDTExtPhDigi> outExtP2Ph;
0882 vector<L1Phase2MuDTThDigi> outP2Th;
0883 vector<L1Phase2MuDTExtThDigi> outExtP2Th;
0884
0885
0886 if (!skip_processing_)
0887 for (auto& ch_correlatedMetaPrimitives : coMetaPrimitives) {
0888 assignIndex(ch_correlatedMetaPrimitives.second);
0889 }
0890
0891 for (auto& ch_correlatedMetaPrimitives : coMetaPrimitives) {
0892 for (const auto& metaPrimitiveIt : ch_correlatedMetaPrimitives.second) {
0893 DTChamberId chId(metaPrimitiveIt.rawId);
0894 DTSuperLayerId slId(metaPrimitiveIt.rawId);
0895 if (debug_)
0896 LogDebug("DTTrigPhase2Prod") << "looping in final vector: SuperLayerId" << chId << " x=" << metaPrimitiveIt.x
0897 << " quality=" << metaPrimitiveIt.quality
0898 << " BX=" << round(metaPrimitiveIt.t0 / 25.) << " index=" << metaPrimitiveIt.index;
0899
0900 int sectorTP = chId.sector();
0901
0902
0903 if (sectorTP == 13)
0904 sectorTP = 4;
0905 if (sectorTP == 14)
0906 sectorTP = 10;
0907 sectorTP = sectorTP - 1;
0908 int sl = 0;
0909 if (metaPrimitiveIt.quality < LOWLOWQ || metaPrimitiveIt.quality == CHIGHQ) {
0910 if (inner(metaPrimitiveIt))
0911 sl = 1;
0912 else
0913 sl = 3;
0914 }
0915
0916 float tp_t0 =
0917 (metaPrimitiveIt.t0 - shift_back * LHC_CLK_FREQ) * ((float)TIME_TO_TDC_COUNTS / (float)LHC_CLK_FREQ);
0918
0919 if (debug_)
0920 LogDebug("DTTrigPhase2Prod") << "pushing back phase-2 dataformat carlo-federica dataformat";
0921
0922 if (slId.superLayer() != 2) {
0923 if (df_extended_ == 1 || df_extended_ == 2) {
0924 int pathWireId[8] = {metaPrimitiveIt.wi1,
0925 metaPrimitiveIt.wi2,
0926 metaPrimitiveIt.wi3,
0927 metaPrimitiveIt.wi4,
0928 metaPrimitiveIt.wi5,
0929 metaPrimitiveIt.wi6,
0930 metaPrimitiveIt.wi7,
0931 metaPrimitiveIt.wi8};
0932
0933 int pathTDC[8] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
0934 max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
0935 max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
0936 max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1),
0937 max((int)round(metaPrimitiveIt.tdc5 - shift_back * LHC_CLK_FREQ), -1),
0938 max((int)round(metaPrimitiveIt.tdc6 - shift_back * LHC_CLK_FREQ), -1),
0939 max((int)round(metaPrimitiveIt.tdc7 - shift_back * LHC_CLK_FREQ), -1),
0940 max((int)round(metaPrimitiveIt.tdc8 - shift_back * LHC_CLK_FREQ), -1)};
0941
0942 int pathLat[8] = {metaPrimitiveIt.lat1,
0943 metaPrimitiveIt.lat2,
0944 metaPrimitiveIt.lat3,
0945 metaPrimitiveIt.lat4,
0946 metaPrimitiveIt.lat5,
0947 metaPrimitiveIt.lat6,
0948 metaPrimitiveIt.lat7,
0949 metaPrimitiveIt.lat8};
0950
0951
0952 outExtP2Ph.emplace_back(
0953 L1Phase2MuDTExtPhDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0954 chId.wheel(),
0955 sectorTP,
0956 chId.station(),
0957 sl,
0958 (int)round(metaPrimitiveIt.phi * PHIRES_CONV),
0959 (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV),
0960 metaPrimitiveIt.quality,
0961 metaPrimitiveIt.index,
0962 tp_t0,
0963 (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),
0964 (int)round(metaPrimitiveIt.x * 1000),
0965 (int)round(metaPrimitiveIt.tanPhi * 1000),
0966 (int)round(metaPrimitiveIt.phi_cmssw * PHIRES_CONV),
0967 (int)round(metaPrimitiveIt.phiB_cmssw * PHIBRES_CONV),
0968 metaPrimitiveIt.rpcFlag,
0969 pathWireId,
0970 pathTDC,
0971 pathLat));
0972 }
0973 if (df_extended_ == 0 || df_extended_ == 2) {
0974
0975 outP2Ph.push_back(L1Phase2MuDTPhDigi(
0976 (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
0977 chId.wheel(),
0978 sectorTP,
0979 chId.station(),
0980 sl,
0981 (int)round(metaPrimitiveIt.phi * PHIRES_CONV),
0982 (int)round(metaPrimitiveIt.phiB * PHIBRES_CONV),
0983 metaPrimitiveIt.quality,
0984 metaPrimitiveIt.index,
0985 tp_t0,
0986 (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),
0987 metaPrimitiveIt.rpcFlag
0988 ));
0989 }
0990 } else {
0991 if (df_extended_ == 1 || df_extended_ == 2) {
0992 int pathWireId[4] = {metaPrimitiveIt.wi1, metaPrimitiveIt.wi2, metaPrimitiveIt.wi3, metaPrimitiveIt.wi4};
0993
0994 int pathTDC[4] = {max((int)round(metaPrimitiveIt.tdc1 - shift_back * LHC_CLK_FREQ), -1),
0995 max((int)round(metaPrimitiveIt.tdc2 - shift_back * LHC_CLK_FREQ), -1),
0996 max((int)round(metaPrimitiveIt.tdc3 - shift_back * LHC_CLK_FREQ), -1),
0997 max((int)round(metaPrimitiveIt.tdc4 - shift_back * LHC_CLK_FREQ), -1)};
0998
0999 int pathLat[4] = {metaPrimitiveIt.lat1, metaPrimitiveIt.lat2, metaPrimitiveIt.lat3, metaPrimitiveIt.lat4};
1000
1001
1002 outExtP2Th.emplace_back(
1003 L1Phase2MuDTExtThDigi((int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
1004 chId.wheel(),
1005 sectorTP,
1006 chId.station(),
1007 (int)round(metaPrimitiveIt.phi * ZRES_CONV),
1008 (int)round(metaPrimitiveIt.phiB * KRES_CONV),
1009 metaPrimitiveIt.quality,
1010 metaPrimitiveIt.index,
1011 tp_t0,
1012 (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),
1013 (int)round(metaPrimitiveIt.x * 1000),
1014 (int)round(metaPrimitiveIt.phi_cmssw * ZRES_CONV),
1015 (int)round(metaPrimitiveIt.phiB_cmssw * KRES_CONV),
1016 metaPrimitiveIt.rpcFlag,
1017 pathWireId,
1018 pathTDC,
1019 pathLat));
1020 }
1021 if (df_extended_ == 0 || df_extended_ == 2) {
1022
1023 outP2Th.push_back(L1Phase2MuDTThDigi(
1024 (int)round(metaPrimitiveIt.t0 / (float)LHC_CLK_FREQ) - shift_back,
1025 chId.wheel(),
1026 sectorTP,
1027 chId.station(),
1028 (int)round(metaPrimitiveIt.phi * ZRES_CONV),
1029 (int)round(metaPrimitiveIt.phiB * KRES_CONV),
1030 metaPrimitiveIt.quality,
1031 metaPrimitiveIt.index,
1032 tp_t0,
1033 (int)round(metaPrimitiveIt.chi2 * CHI2RES_CONV),
1034 metaPrimitiveIt.rpcFlag
1035 ));
1036 }
1037 }
1038 }
1039 }
1040
1041
1042 if (useRPC_) {
1043 for (auto rpc_dt_digi = rpc_integrator_->rpcRecHits_translated_.begin();
1044 rpc_dt_digi != rpc_integrator_->rpcRecHits_translated_.end();
1045 rpc_dt_digi++) {
1046 outP2Ph.push_back(*rpc_dt_digi);
1047 }
1048 }
1049
1050
1051 if (df_extended_ == 1 || df_extended_ == 2) {
1052 std::unique_ptr<L1Phase2MuDTExtPhContainer> resultExtP2Ph(new L1Phase2MuDTExtPhContainer);
1053 resultExtP2Ph->setContainer(outExtP2Ph);
1054 iEvent.put(std::move(resultExtP2Ph));
1055 }
1056 if (df_extended_ == 0 || df_extended_ == 2) {
1057 std::unique_ptr<L1Phase2MuDTPhContainer> resultP2Ph(new L1Phase2MuDTPhContainer);
1058 resultP2Ph->setContainer(outP2Ph);
1059 iEvent.put(std::move(resultP2Ph));
1060 }
1061 outExtP2Ph.clear();
1062 outExtP2Ph.erase(outExtP2Ph.begin(), outExtP2Ph.end());
1063 outP2Ph.clear();
1064 outP2Ph.erase(outP2Ph.begin(), outP2Ph.end());
1065
1066
1067 if (df_extended_ == 1 || df_extended_ == 2) {
1068 std::unique_ptr<L1Phase2MuDTExtThContainer> resultExtP2Th(new L1Phase2MuDTExtThContainer);
1069 resultExtP2Th->setContainer(outExtP2Th);
1070 iEvent.put(std::move(resultExtP2Th));
1071 }
1072 if (df_extended_ == 0 || df_extended_ == 2) {
1073 std::unique_ptr<L1Phase2MuDTThContainer> resultP2Th(new L1Phase2MuDTThContainer);
1074 resultP2Th->setContainer(outP2Th);
1075 iEvent.put(std::move(resultP2Th));
1076 }
1077 outExtP2Th.clear();
1078 outExtP2Th.erase(outExtP2Th.begin(), outExtP2Th.end());
1079 outP2Th.clear();
1080 outP2Th.erase(outP2Th.begin(), outP2Th.end());
1081 }
1082
1083 void DTTrigPhase2Prod::endRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
1084 grouping_obj_->finish();
1085 mpathanalyzer_->finish();
1086 mpathqualityenhancer_->finish();
1087 mpathqualityenhancerbayes_->finish();
1088 mpathredundantfilter_->finish();
1089 mpathhitsfilter_->finish();
1090 mpathassociator_->finish();
1091 rpc_integrator_->finish();
1092 };
1093
1094 bool DTTrigPhase2Prod::outer(const metaPrimitive& mp) const {
1095 int counter = (mp.wi5 != -1) + (mp.wi6 != -1) + (mp.wi7 != -1) + (mp.wi8 != -1);
1096 return (counter > 2);
1097 }
1098
1099 bool DTTrigPhase2Prod::inner(const metaPrimitive& mp) const {
1100 int counter = (mp.wi1 != -1) + (mp.wi2 != -1) + (mp.wi3 != -1) + (mp.wi4 != -1);
1101 return (counter > 2);
1102 }
1103
1104 bool DTTrigPhase2Prod::hasPosRF(int wh, int sec) const { return wh > 0 || (wh == 0 && sec % 4 > 1); }
1105
1106 void DTTrigPhase2Prod::printmP(const string& ss, const metaPrimitive& mP) const {
1107 DTSuperLayerId slId(mP.rawId);
1108 LogInfo("DTTrigPhase2Prod") << ss << (int)slId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
1109 << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
1110 << setw(5) << left << mP.tdc1 << " " << setw(5) << left << mP.tdc2 << " " << setw(5)
1111 << left << mP.tdc3 << " " << setw(5) << left << mP.tdc4 << " " << setw(10) << right
1112 << mP.x << " " << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " "
1113 << setw(13) << left << mP.chi2 << " r:" << rango(mP);
1114 }
1115
1116 void DTTrigPhase2Prod::printmP(const metaPrimitive& mP) const {
1117 DTSuperLayerId slId(mP.rawId);
1118 LogInfo("DTTrigPhase2Prod") << (int)slId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left << mP.wi2
1119 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " " << setw(5)
1120 << left << mP.tdc1 << " " << setw(5) << left << mP.tdc2 << " " << setw(5) << left
1121 << mP.tdc3 << " " << setw(5) << left << mP.tdc4 << " " << setw(10) << right << mP.x << " "
1122 << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13)
1123 << left << mP.chi2 << " r:" << rango(mP) << std::endl;
1124 }
1125
1126 void DTTrigPhase2Prod::printmPC(const string& ss, const metaPrimitive& mP) const {
1127 DTChamberId ChId(mP.rawId);
1128 LogInfo("DTTrigPhase2Prod") << ss << (int)ChId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left
1129 << mP.wi2 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " "
1130 << setw(2) << left << mP.wi5 << " " << setw(2) << left << mP.wi6 << " " << setw(2) << left
1131 << mP.wi7 << " " << setw(2) << left << mP.wi8 << " " << setw(5) << left << mP.tdc1 << " "
1132 << setw(5) << left << mP.tdc2 << " " << setw(5) << left << mP.tdc3 << " " << setw(5)
1133 << left << mP.tdc4 << " " << setw(5) << left << mP.tdc5 << " " << setw(5) << left
1134 << mP.tdc6 << " " << setw(5) << left << mP.tdc7 << " " << setw(5) << left << mP.tdc8
1135 << " " << setw(2) << left << mP.lat1 << " " << setw(2) << left << mP.lat2 << " "
1136 << setw(2) << left << mP.lat3 << " " << setw(2) << left << mP.lat4 << " " << setw(2)
1137 << left << mP.lat5 << " " << setw(2) << left << mP.lat6 << " " << setw(2) << left
1138 << mP.lat7 << " " << setw(2) << left << mP.lat8 << " " << setw(10) << right << mP.x << " "
1139 << setw(9) << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13)
1140 << left << mP.chi2 << " r:" << rango(mP);
1141 }
1142
1143 void DTTrigPhase2Prod::printmPC(const metaPrimitive& mP) const {
1144 DTChamberId ChId(mP.rawId);
1145 LogInfo("DTTrigPhase2Prod") << (int)ChId << "\t " << setw(2) << left << mP.wi1 << " " << setw(2) << left << mP.wi2
1146 << " " << setw(2) << left << mP.wi3 << " " << setw(2) << left << mP.wi4 << " " << setw(2)
1147 << left << mP.wi5 << " " << setw(2) << left << mP.wi6 << " " << setw(2) << left << mP.wi7
1148 << " " << setw(2) << left << mP.wi8 << " " << setw(5) << left << mP.tdc1 << " " << setw(5)
1149 << left << mP.tdc2 << " " << setw(5) << left << mP.tdc3 << " " << setw(5) << left
1150 << mP.tdc4 << " " << setw(5) << left << mP.tdc5 << " " << setw(5) << left << mP.tdc6
1151 << " " << setw(5) << left << mP.tdc7 << " " << setw(5) << left << mP.tdc8 << " "
1152 << setw(2) << left << mP.lat1 << " " << setw(2) << left << mP.lat2 << " " << setw(2)
1153 << left << mP.lat3 << " " << setw(2) << left << mP.lat4 << " " << setw(2) << left
1154 << mP.lat5 << " " << setw(2) << left << mP.lat6 << " " << setw(2) << left << mP.lat7
1155 << " " << setw(2) << left << mP.lat8 << " " << setw(10) << right << mP.x << " " << setw(9)
1156 << left << mP.tanPhi << " " << setw(5) << left << mP.t0 << " " << setw(13) << left
1157 << mP.chi2 << " r:" << rango(mP) << std::endl;
1158 }
1159
1160 int DTTrigPhase2Prod::rango(const metaPrimitive& mp) const {
1161 if (mp.quality == 1 or mp.quality == 2)
1162 return 3;
1163 if (mp.quality == 3 or mp.quality == 4)
1164 return 4;
1165 return mp.quality;
1166 }
1167
1168 void DTTrigPhase2Prod::assignIndex(std::vector<metaPrimitive>& inMPaths) {
1169 std::map<int, std::vector<metaPrimitive>> primsPerBX;
1170 for (const auto& metaPrimitive : inMPaths) {
1171 int BX = round(metaPrimitive.t0 / 25.);
1172 primsPerBX[BX].push_back(metaPrimitive);
1173 }
1174 inMPaths.clear();
1175 for (auto& prims : primsPerBX) {
1176 assignIndexPerBX(prims.second);
1177 for (const auto& primitive : prims.second)
1178 if (primitive.index <= max_index_)
1179 inMPaths.push_back(primitive);
1180 }
1181 }
1182
1183 void DTTrigPhase2Prod::assignIndexPerBX(std::vector<metaPrimitive>& inMPaths) {
1184
1185 uint32_t rawId = -1;
1186 int numP = -1;
1187 for (auto& metaPrimitiveIt : inMPaths) {
1188 numP++;
1189 rawId = metaPrimitiveIt.rawId;
1190 int iOrder = assignQualityOrder(metaPrimitiveIt);
1191 int inf = 0;
1192 int numP2 = -1;
1193 for (auto& metaPrimitiveItN : inMPaths) {
1194 int nOrder = assignQualityOrder(metaPrimitiveItN);
1195 numP2++;
1196 if (rawId != metaPrimitiveItN.rawId)
1197 continue;
1198 if (numP2 == numP) {
1199 metaPrimitiveIt.index = inf;
1200 break;
1201 } else if (iOrder < nOrder) {
1202 inf++;
1203 } else if (iOrder > nOrder) {
1204 metaPrimitiveItN.index++;
1205 } else if (iOrder == nOrder) {
1206 if (std::abs(metaPrimitiveIt.phiB) >= std::abs(metaPrimitiveItN.phiB)) {
1207 inf++;
1208 } else if (std::abs(metaPrimitiveIt.phiB) < std::abs(metaPrimitiveItN.phiB)) {
1209 metaPrimitiveItN.index++;
1210 }
1211 }
1212 }
1213 }
1214 }
1215
1216 int DTTrigPhase2Prod::assignQualityOrder(const metaPrimitive& mP) const {
1217 if (mP.quality > 8 || mP.quality < 1)
1218 return -1;
1219
1220 return qmap_.find(mP.quality)->second;
1221 }
1222
1223 std::vector<DTDigiCollection*> DTTrigPhase2Prod::distribDigis(std::queue<std::pair<DTLayerId, DTDigi>>& inQ) {
1224 std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*> tmpVector;
1225 tmpVector.clear();
1226 std::vector<DTDigiCollection*> collVector;
1227 collVector.clear();
1228 while (!inQ.empty()) {
1229 processDigi(inQ, tmpVector);
1230 }
1231 for (auto& sQ : tmpVector) {
1232 DTDigiCollection tmpColl;
1233 while (!sQ->empty()) {
1234 tmpColl.insertDigi((sQ->front().first), (sQ->front().second));
1235 sQ->pop();
1236 }
1237 collVector.push_back(&tmpColl);
1238 }
1239 return collVector;
1240 }
1241
1242 void DTTrigPhase2Prod::processDigi(std::queue<std::pair<DTLayerId, DTDigi>>& inQ,
1243 std::vector<std::queue<std::pair<DTLayerId, DTDigi>>*>& vec) {
1244 bool classified = false;
1245 if (!vec.empty()) {
1246 for (auto& sC : vec) {
1247 if ((sC->front().second.time() + superCelltimewidth_) > inQ.front().second.time()) {
1248
1249 if (TMath::Abs(sC->front().second.wire() - inQ.front().second.wire()) <= superCellhalfspacewidth_) {
1250
1251 sC->push(std::move(inQ.front()));
1252 classified = true;
1253 }
1254 }
1255 }
1256 }
1257 if (classified) {
1258 inQ.pop();
1259 return;
1260 }
1261
1262 std::queue<std::pair<DTLayerId, DTDigi>> newQueue;
1263
1264 std::pair<DTLayerId, DTDigi> tmpPair;
1265 tmpPair = std::move(inQ.front());
1266 newQueue.push(tmpPair);
1267 inQ.pop();
1268
1269 vec.push_back(&newQueue);
1270 }
1271
1272 void DTTrigPhase2Prod::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
1273
1274 edm::ParameterSetDescription desc;
1275 desc.add<edm::InputTag>("digiTag", edm::InputTag("CalibratedDigis"));
1276 desc.add<int>("timeTolerance", 999999);
1277 desc.add<double>("tanPhiTh", 1.0);
1278 desc.add<double>("tanPhiThw2max", 1.3);
1279 desc.add<double>("tanPhiThw2min", 0.5);
1280 desc.add<double>("tanPhiThw1max", 0.9);
1281 desc.add<double>("tanPhiThw1min", 0.2);
1282 desc.add<double>("tanPhiThw0", 0.5);
1283 desc.add<double>("chi2Th", 0.01);
1284 desc.add<double>("chi2corTh", 0.1);
1285 desc.add<bool>("useBX_correlation", false);
1286 desc.add<double>("dT0_correlate_TP", 25.0);
1287 desc.add<int>("dBX_correlate_TP", 0);
1288 desc.add<double>("dTanPsi_correlate_TP", 99999.0);
1289 desc.add<bool>("clean_chi2_correlation", true);
1290 desc.add<bool>("allow_confirmation", true);
1291 desc.add<double>("minx_match_2digis", 1.0);
1292 desc.add<int>("scenario", 0);
1293 desc.add<int>("df_extended", 0);
1294 desc.add<int>("co_option", 0);
1295 desc.add<int>("co_quality", 0);
1296 desc.add<int>("max_primitives", 999);
1297 desc.add<bool>("output_mixer", false);
1298 desc.add<bool>("output_latpredictor", false);
1299 desc.add<bool>("output_slfitter", false);
1300 desc.add<bool>("output_slfilter", false);
1301 desc.add<bool>("output_confirmed", false);
1302 desc.add<bool>("output_matcher", false);
1303 desc.add<edm::FileInPath>("ttrig_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_ttrig.txt"));
1304 desc.add<edm::FileInPath>("z_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_z.txt"));
1305 desc.add<edm::FileInPath>("lut_sl1", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_sl1.dat"));
1306 desc.add<edm::FileInPath>("lut_sl2", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_slx.dat"));
1307 desc.add<edm::FileInPath>("lut_sl3", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_sl3.dat"));
1308 desc.add<edm::FileInPath>("lut_2sl", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/fitterlut_2sl.dat"));
1309 desc.add<edm::FileInPath>("shift_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/wire_rawId_x.txt"));
1310 desc.add<edm::FileInPath>("maxdrift_filename",
1311 edm::FileInPath("L1Trigger/DTTriggerPhase2/data/drift_time_per_chamber.txt"));
1312 desc.add<edm::FileInPath>("shift_theta_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/theta_shift.txt"));
1313 desc.add<edm::FileInPath>("global_coords_filename",
1314 edm::FileInPath("L1Trigger/DTTriggerPhase2/data/global_coord_perp_x_phi0.txt"));
1315 desc.add<edm::FileInPath>("laterality_filename",
1316 edm::FileInPath("L1Trigger/DTTriggerPhase2/data/lat_predictions.dat"));
1317 desc.add<int>("algo", 0);
1318 desc.add<int>("minHits4Fit", 3);
1319 desc.add<bool>("splitPathPerSL", true);
1320 desc.addUntracked<bool>("debug", false);
1321 desc.addUntracked<bool>("dump", false);
1322 desc.add<edm::InputTag>("rpcRecHits", edm::InputTag("rpcRecHits"));
1323 desc.add<bool>("useRPC", false);
1324 desc.add<int>("bx_window", 1);
1325 desc.add<double>("phi_window", 50.0);
1326 desc.add<int>("max_quality_to_overwrite_t0", 9);
1327 desc.add<bool>("storeAllRPCHits", false);
1328 desc.add<bool>("activateBuffer", false);
1329 desc.add<double>("superCelltimewidth", 400);
1330 desc.add<int>("superCellspacewidth", 20);
1331 {
1332 edm::ParameterSetDescription psd0;
1333 psd0.addUntracked<bool>("debug", false);
1334 psd0.add<double>("angletan", 0.3);
1335 psd0.add<double>("anglebinwidth", 1.0);
1336 psd0.add<double>("posbinwidth", 2.1);
1337 psd0.add<double>("maxdeltaAngDeg", 10);
1338 psd0.add<double>("maxdeltaPos", 10);
1339 psd0.add<int>("UpperNumber", 6);
1340 psd0.add<int>("LowerNumber", 4);
1341 psd0.add<double>("MaxDistanceToWire", 0.03);
1342 psd0.add<int>("minNLayerHits", 6);
1343 psd0.add<int>("minSingleSLHitsMax", 3);
1344 psd0.add<int>("minSingleSLHitsMin", 3);
1345 psd0.add<bool>("allowUncorrelatedPatterns", true);
1346 psd0.add<int>("minUncorrelatedHits", 3);
1347 desc.add<edm::ParameterSetDescription>("HoughGrouping", psd0);
1348 }
1349 {
1350 edm::ParameterSetDescription psd0;
1351 psd0.add<edm::FileInPath>(
1352 "pattern_filename", edm::FileInPath("L1Trigger/DTTriggerPhase2/data/PseudoBayesPatterns_uncorrelated_v0.root"));
1353 psd0.addUntracked<bool>("debug", false);
1354 psd0.add<int>("minNLayerHits", 3);
1355 psd0.add<int>("minSingleSLHitsMax", 3);
1356 psd0.add<int>("minSingleSLHitsMin", 0);
1357 psd0.add<int>("allowedVariance", 1);
1358 psd0.add<bool>("allowDuplicates", false);
1359 psd0.add<bool>("setLateralities", true);
1360 psd0.add<bool>("allowUncorrelatedPatterns", true);
1361 psd0.add<int>("minUncorrelatedHits", 3);
1362 psd0.add<bool>("saveOnPlace", true);
1363 psd0.add<int>("maxPathsPerMatch", 256);
1364 desc.add<edm::ParameterSetDescription>("PseudoBayesPattern", psd0);
1365 }
1366 descriptions.add("dtTriggerPhase2PrimitiveDigis", desc);
1367 }
1368
1369 DEFINE_FWK_MODULE(DTTrigPhase2Prod);