Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:39

0001 // $Id:
0002 
0003 //-----------------------------------------------------------------------------
0004 // Implementation file for class : RPCTechnicalTrigger
0005 //
0006 // 2008-10-15 : Andres Osorio
0007 //-----------------------------------------------------------------------------
0008 
0009 //=============================================================================
0010 // Standard constructor, initializes variables
0011 //=============================================================================
0012 
0013 // Include files
0014 
0015 // local
0016 #include "L1Trigger/RPCTechnicalTrigger/interface/RPCTechnicalTrigger.h"
0017 #include "L1Trigger/RPCTechnicalTrigger/interface/ProcessTestSignal.h"
0018 #include "L1Trigger/RPCTechnicalTrigger/interface/RBCProcessRPCDigis.h"
0019 #include "L1Trigger/RPCTechnicalTrigger/interface/RBCProcessRPCSimDigis.h"
0020 
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/MessageLogger/interface/MessageDrop.h"
0023 
0024 //=============================================================================
0025 // Standard constructor, initializes variables
0026 //=============================================================================
0027 
0028 namespace {
0029   //...........................................................................
0030   //For the pointing Logic: declare here the first sector of each quadrant
0031   //
0032   constexpr std::array<int, 10> s_quadrants = {{2, 3, 4, 5, 6, 7, 8, 9, 10, 11}};
0033 
0034   //The wheelTtu is addressed using -2, -1, 0, 1, 2
0035   constexpr unsigned int kWheelOffset = 2;
0036   constexpr std::array<int, 5> wheelTtu = {{3, 3, 2, 1, 1}};
0037 }  // namespace
0038 
0039 RPCTechnicalTrigger::RPCTechnicalTrigger(const edm::ParameterSet& iConfig)
0040     : m_verbosity{iConfig.getUntrackedParameter<int>("Verbosity", 0)},
0041       m_useEventSetup{iConfig.getUntrackedParameter<int>("UseEventSetup", 0)},
0042       m_ttBits{iConfig.getParameter<std::vector<unsigned>>("BitNumbers")},
0043       m_ttNames{iConfig.getParameter<std::vector<std::string>>("BitNames")},
0044       m_rpcDigiLabel{iConfig.getParameter<edm::InputTag>("RPCDigiLabel")},
0045       m_rpcDigiToken{consumes<RPCDigiCollection>(m_rpcDigiLabel)},
0046       m_useRPCSimLink{iConfig.getUntrackedParameter<int>("UseRPCSimLink", 0)},
0047       m_rpcGeometryToken(esConsumes<RPCGeometry, MuonGeometryRecord>()) {
0048   std::string configFile = iConfig.getParameter<std::string>("ConfigFile");
0049 
0050   edm::FileInPath f1("L1Trigger/RPCTechnicalTrigger/data/" + configFile);
0051   m_configFile = f1.fullPath();
0052 
0053   if (m_verbosity) {
0054     LogTrace("RPCTechnicalTrigger") << m_rpcDigiLabel << '\n' << std::endl;
0055 
0056     LogTrace("RPCTechnicalTrigger") << "\nConfiguration file used for UseEventSetup = 0 \n"
0057                                     << m_configFile << '\n'
0058                                     << std::endl;
0059   }
0060 
0061   //...........................................................................
0062   //... There are three Technical Trigger Units Boards: 1 can handle 2 Wheels
0063   //... n_Wheels sets the number of wheels attached to board with index boardIndex
0064 
0065   constexpr std::array<int, 3> boardIndex = {{1, 2, 3}};
0066   constexpr std::array<int, 3> nWheels = {{2, 1, 2}};
0067 
0068   m_ttu[0] = TTUEmulator(boardIndex[0], nWheels[0]);
0069   m_ttu[1] = TTUEmulator(boardIndex[1], nWheels[1]);
0070   m_ttu[2] = TTUEmulator(boardIndex[2], nWheels[2]);
0071 
0072   //... This is second line that delivers in parallel a second trigger
0073   m_ttuRbcLine[0] = TTUEmulator(boardIndex[0], nWheels[0]);
0074   m_ttuRbcLine[1] = TTUEmulator(boardIndex[1], nWheels[1]);
0075   m_ttuRbcLine[2] = TTUEmulator(boardIndex[2], nWheels[2]);
0076 
0077   //...........................................................................
0078 
0079   m_hasConfig = false;
0080   produces<L1GtTechnicalTriggerRecord>();
0081   consumes<edm::DetSetVector<RPCDigiSimLink>>(edm::InputTag("simMuonRPCDigis", "RPCDigiSimLink", ""));
0082   if (m_useEventSetup >= 1) {
0083     m_pRBCSpecsToken = esConsumes<RBCBoardSpecs, RBCBoardSpecsRcd, edm::Transition::BeginRun>();
0084     m_pTTUSpecsToken = esConsumes<TTUBoardSpecs, TTUBoardSpecsRcd, edm::Transition::BeginRun>();
0085   }
0086 }
0087 
0088 RPCTechnicalTrigger::~RPCTechnicalTrigger() {}
0089 
0090 //=============================================================================
0091 void RPCTechnicalTrigger::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0092   bool status(false);
0093 
0094   edm::Handle<RPCDigiCollection> pIn;
0095 
0096   edm::Handle<edm::DetSetVector<RPCDigiSimLink>> simIn;
0097 
0098   std::unique_ptr<L1GtTechnicalTriggerRecord> output(new L1GtTechnicalTriggerRecord());
0099 
0100   //.   Set up RPC geometry
0101   edm::ESHandle<RPCGeometry> rpcGeometry = iSetup.getHandle(m_rpcGeometryToken);
0102 
0103   std::unique_ptr<ProcessInputSignal> signal;
0104   if (m_useRPCSimLink == 0) {
0105     iEvent.getByToken(m_rpcDigiToken, pIn);
0106     if (!pIn.isValid()) {
0107       edm::LogError("RPCTechnicalTrigger") << "can't find RPCDigiCollection with label: " << m_rpcDigiLabel << '\n';
0108       iEvent.put(std::move(output));
0109       return;
0110     }
0111 
0112     signal = std::make_unique<RBCProcessRPCDigis>(rpcGeometry, pIn);
0113 
0114   } else {
0115     iEvent.getByLabel("simMuonRPCDigis", "RPCDigiSimLink", simIn);
0116 
0117     if (!simIn.isValid()) {
0118       edm::LogError("RPCTechnicalTrigger") << "can't find RPCDigiCollection with label: " << m_rpcDigiLabel << '\n';
0119       iEvent.put(std::move(output));
0120       return;
0121     }
0122     signal = std::make_unique<RBCProcessRPCSimDigis>(rpcGeometry, simIn);
0123   }
0124 
0125   LogDebug("RPCTechnicalTrigger") << "signal object created" << '\n';
0126 
0127   if (!m_hasConfig) {
0128     edm::LogError("RPCTechnicalTrigger") << "cannot read hardware configuration \n";
0129     iEvent.put(std::move(output));
0130     return;
0131   }
0132 
0133   status = signal->next();
0134 
0135   if (!status) {
0136     iEvent.put(std::move(output));
0137     return;
0138   }
0139 
0140   auto* input = signal->retrievedata();
0141 
0142   std::vector<L1GtTechnicalTrigger> ttVec(m_ttBits.size());
0143 
0144   //. distribute data to different TTU emulator instances and process it
0145   std::bitset<5> triggerbits;
0146 
0147   std::vector<std::unique_ptr<TTUResults>> serializedInfoLine1;
0148   std::vector<std::unique_ptr<TTUResults>> serializedInfoLine2;
0149 
0150   for (int k = 0; k < kMaxTtuBoards; ++k) {
0151     m_ttu[k].processTtu(input);
0152 
0153     //work out Pointing Logic to Tracker
0154     for (auto firstSector : s_quadrants)
0155       m_ttuRbcLine[k].processTtu(input, firstSector);
0156 
0157     //...for trigger 1
0158     for (auto const& out : m_ttu[k].m_triggerBxVec)
0159       serializedInfoLine1.emplace_back(std::make_unique<TTUResults>(k, out.m_bx, out.m_trigger[0], out.m_trigger[1]));
0160     m_ttu[k].clearTriggerResponse();
0161 
0162     //...for trigger 2
0163     for (auto const& out : m_ttuRbcLine[k].m_triggerBxVec)
0164       serializedInfoLine2.push_back(
0165           std::make_unique<TTUResults>(k, out.m_bx, out.m_trigger[0], out.m_trigger[1], out.m_wedge));
0166 
0167     m_ttuRbcLine[k].clearTriggerResponse();
0168   }
0169 
0170   //.. write results to technical trigger bits
0171   int bx(0);
0172   int infoSize(0);
0173 
0174   infoSize = serializedInfoLine1.size();
0175 
0176   auto sortByBx = [](auto& iLHS, auto& iRHS) { return iLHS->m_bx < iRHS->m_bx; };
0177   std::sort(serializedInfoLine1.begin(), serializedInfoLine1.end(), sortByBx);
0178 
0179   if (m_verbosity) {
0180     for (auto& ttu : serializedInfoLine1) {
0181       if (abs(ttu->m_bx) <= 1)
0182         std::cout << "RPCTechnicalTrigger> " << ttu->m_ttuidx << '\t' << ttu->m_bx << '\t' << ttu->m_trigWheel1 << '\t'
0183                   << ttu->m_trigWheel2 << '\n';
0184     }
0185   }
0186 
0187   bool has_bx0 = false;
0188 
0189   for (int k = 0; k < infoSize; k += kMaxTtuBoards) {
0190     bx = serializedInfoLine1[k]->m_bx;
0191 
0192     if (bx == 0) {
0193       triggerbits.set(0, serializedInfoLine1[k]->m_trigWheel2);
0194       triggerbits.set(1, serializedInfoLine1[k]->m_trigWheel1);
0195       triggerbits.set(2, serializedInfoLine1[k + 1]->m_trigWheel1);
0196       triggerbits.set(3, serializedInfoLine1[k + 2]->m_trigWheel1);
0197       triggerbits.set(4, serializedInfoLine1[k + 2]->m_trigWheel2);
0198 
0199       bool five_wheels_OR = triggerbits.any();
0200 
0201       ttVec.at(0) = L1GtTechnicalTrigger(
0202           m_ttNames.at(0), m_ttBits.at(0), bx, five_wheels_OR);  // bit 24 = Or 5 wheels in TTU mode
0203       ttVec.at(2) = L1GtTechnicalTrigger(m_ttNames.at(2), m_ttBits.at(2), bx, triggerbits[0]);  // bit 26
0204       ttVec.at(3) = L1GtTechnicalTrigger(m_ttNames.at(3), m_ttBits.at(3), bx, triggerbits[1]);  // bit 27
0205       ttVec.at(4) = L1GtTechnicalTrigger(m_ttNames.at(4), m_ttBits.at(4), bx, triggerbits[2]);  // bit 28
0206       ttVec.at(5) = L1GtTechnicalTrigger(m_ttNames.at(5), m_ttBits.at(5), bx, triggerbits[3]);  // bit 29
0207       ttVec.at(6) = L1GtTechnicalTrigger(m_ttNames.at(6), m_ttBits.at(6), bx, triggerbits[4]);  // bit 30
0208 
0209       triggerbits.reset();
0210 
0211       has_bx0 = true;
0212 
0213       break;
0214 
0215     } else
0216       continue;
0217   }
0218 
0219   infoSize = serializedInfoLine2.size();
0220 
0221   std::sort(serializedInfoLine2.begin(), serializedInfoLine2.end(), sortByBx);
0222 
0223   if (m_verbosity) {
0224     for (auto& ttu : serializedInfoLine2) {
0225       if (abs(ttu->m_bx) <= 1)
0226         std::cout << "RPCTechnicalTrigger> " << ttu->m_ttuidx << '\t' << ttu->m_bx << '\t' << ttu->m_trigWheel1 << '\t'
0227                   << ttu->m_trigWheel2 << '\t' << ttu->m_wedge << '\n';
0228     }
0229   }
0230 
0231   auto ttuResultsByQuadrant = convertToMap(serializedInfoLine2);
0232 
0233   std::bitset<8> triggerCoincidence;
0234   triggerCoincidence.reset();
0235 
0236   // searchCoincidence( W-2 , W0 )
0237   bool result = searchCoincidence(-2, 0, ttuResultsByQuadrant);
0238   triggerCoincidence.set(0, result);
0239 
0240   // searchCoincidence( W-2 , W+1 )
0241   result = searchCoincidence(-2, 1, ttuResultsByQuadrant);
0242   triggerCoincidence.set(1, result);
0243 
0244   // searchCoincidence( W-1 , W0  )
0245   result = searchCoincidence(-1, 0, ttuResultsByQuadrant);
0246   triggerCoincidence.set(2, result);
0247 
0248   // searchCoincidence( W-1 , W+1 )
0249   result = searchCoincidence(-1, 1, ttuResultsByQuadrant);
0250   triggerCoincidence.set(3, result);
0251 
0252   // searchCoincidence( W-1 , W+2 )
0253   result = searchCoincidence(-1, 2, ttuResultsByQuadrant);
0254   triggerCoincidence.set(4, result);
0255 
0256   // searchCoincidence( W0  , W0  )
0257   result = searchCoincidence(0, 0, ttuResultsByQuadrant);
0258   triggerCoincidence.set(5, result);
0259 
0260   // searchCoincidence( W+1 , W0  )
0261   result = searchCoincidence(1, 0, ttuResultsByQuadrant);
0262   triggerCoincidence.set(6, result);
0263 
0264   // searchCoincidence( W+2 , W0  )
0265   result = searchCoincidence(2, 0, ttuResultsByQuadrant);
0266   triggerCoincidence.set(7, result);
0267 
0268   bool five_wheels_OR = triggerCoincidence.any();
0269 
0270   if (m_verbosity)
0271     std::cout << "RPCTechnicalTrigger> pointing trigger: " << five_wheels_OR << '\n';
0272 
0273   ttVec.at(1) =
0274       L1GtTechnicalTrigger(m_ttNames.at(1), m_ttBits.at(1), bx, five_wheels_OR);  // bit 25 = Or 5 wheels in RBC mode
0275 
0276   triggerCoincidence.reset();
0277 
0278   //...check that data appeared at bx=0
0279 
0280   if (!has_bx0) {
0281     iEvent.put(std::move(output));
0282     LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger> end of event loop" << std::endl;
0283     return;
0284   }
0285 
0286   output->setGtTechnicalTrigger(ttVec);
0287   iEvent.put(std::move(output));
0288 
0289   //.... all done
0290 
0291   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger> end of event loop" << std::endl;
0292 }
0293 // ------------ method called once each job just before starting event loop  ------------
0294 void RPCTechnicalTrigger::beginRun(edm::Run const& iRun, const edm::EventSetup& evtSetup) {
0295   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger::beginRun> starts" << std::endl;
0296 
0297   //..  Get Board Specifications (hardware configuration)
0298 
0299   if (m_useEventSetup >= 1) {
0300     edm::ESHandle<RBCBoardSpecs> pRBCSpecs = evtSetup.getHandle(m_pRBCSpecsToken);
0301 
0302     edm::ESHandle<TTUBoardSpecs> pTTUSpecs = evtSetup.getHandle(m_pTTUSpecsToken);
0303 
0304     if (!pRBCSpecs.isValid() || !pTTUSpecs.isValid()) {
0305       edm::LogError("RPCTechnicalTrigger") << "can't find RBC/TTU BoardSpecsRcd" << '\n';
0306       m_hasConfig = false;
0307     } else {
0308       m_rbcspecs = pRBCSpecs.product();
0309       m_ttuspecs = pTTUSpecs.product();
0310       m_hasConfig = true;
0311     }
0312 
0313   } else {
0314     // read hardware configuration from file
0315     m_readConfig = std::make_unique<TTUConfigurator>(m_configFile);
0316 
0317     if (m_readConfig->m_hasConfig) {
0318       m_readConfig->process();
0319       m_rbcspecs = m_readConfig->getRbcSpecs();
0320       m_ttuspecs = m_readConfig->getTtuSpecs();
0321       m_hasConfig = true;
0322     }
0323 
0324     else
0325       m_hasConfig = false;
0326   }
0327 
0328   if (m_hasConfig) {
0329     //... Initialize all
0330 
0331     for (int k = 0; k < kMaxTtuBoards; ++k) {
0332       m_ttu[k].SetLineId(1);
0333       m_ttuRbcLine[k].SetLineId(2);
0334 
0335       m_ttu[k].setSpecifications(m_ttuspecs, m_rbcspecs);
0336       m_ttuRbcLine[k].setSpecifications(m_ttuspecs, m_rbcspecs);
0337 
0338       m_ttu[k].initialise();
0339       m_ttuRbcLine[k].initialise();
0340     }
0341   }
0342 }
0343 
0344 //
0345 std::map<int, RPCTechnicalTrigger::TTUResults*> RPCTechnicalTrigger::convertToMap(
0346     const std::vector<std::unique_ptr<TTUResults>>& ttuResults) const {
0347   std::map<int, TTUResults*> returnValue;
0348   auto itr = ttuResults.begin();
0349 
0350   while (itr != ttuResults.end()) {
0351     if ((*itr)->m_bx != 0) {
0352       ++itr;
0353       continue;
0354     }
0355 
0356     int key(0);
0357     key = 1000 * ((*itr)->m_ttuidx + 1) + 1 * (*itr)->m_wedge;
0358     returnValue[key] = itr->get();
0359     ++itr;
0360   }
0361 
0362   return returnValue;
0363 }
0364 
0365 //...RBC pointing logic to tracker bit 25: hardwired
0366 bool RPCTechnicalTrigger::searchCoincidence(int wheel1,
0367                                             int wheel2,
0368                                             std::map<int, TTUResults*> const& ttuResultsByQuadrant) const {
0369   std::map<int, TTUResults*>::const_iterator itr;
0370   bool topRight(false);
0371   bool botLeft(false);
0372 
0373   int indxW1 = wheelTtu[wheel1 + kWheelOffset];
0374   int indxW2 = wheelTtu[wheel2 + kWheelOffset];
0375 
0376   int k(0);
0377   int key(0);
0378   bool finalTrigger(false);
0379   int maxTopQuadrants = 4;
0380 
0381   //work out Pointing Logic to Tracker
0382 
0383   for (auto firstSector : s_quadrants) {
0384     key = 1000 * (indxW1) + firstSector;
0385 
0386     itr = ttuResultsByQuadrant.find(key);
0387     if (itr != ttuResultsByQuadrant.end())
0388       topRight = (*itr).second->getTriggerForWheel(wheel1);
0389 
0390     //std::cout << "W1: " << wheel1 << " " << "sec: " << firstSector << " dec: " << topRight << '\n';
0391 
0392     key = 1000 * (indxW2) + firstSector + 5;
0393 
0394     itr = ttuResultsByQuadrant.find(key);
0395 
0396     if (itr != ttuResultsByQuadrant.end())
0397       botLeft = (*itr).second->getTriggerForWheel(wheel2);
0398 
0399     //std::cout << "W2: " << wheel2 << " " << "sec: " << firstSector + 5 << " dec: " << botLeft << '\n';
0400 
0401     finalTrigger |= (topRight && botLeft);
0402 
0403     ++k;
0404 
0405     if (k > maxTopQuadrants)
0406       break;
0407   }
0408 
0409   //Try the opposite now
0410 
0411   k = 0;
0412 
0413   for (auto firstSector : s_quadrants) {
0414     key = 1000 * (indxW2) + firstSector;
0415 
0416     itr = ttuResultsByQuadrant.find(key);
0417     if (itr != ttuResultsByQuadrant.end())
0418       topRight = (*itr).second->getTriggerForWheel(wheel1);
0419 
0420     //std::cout << "W1: " << wheel1 << " " << "sec: " << firstSector << " dec: " << topRight << '\n';
0421 
0422     key = 1000 * (indxW1) + firstSector + 5;
0423 
0424     itr = ttuResultsByQuadrant.find(key);
0425 
0426     if (itr != ttuResultsByQuadrant.end())
0427       botLeft = (*itr).second->getTriggerForWheel(wheel2);
0428 
0429     //std::cout << "W2: " << wheel2 << " " << "sec: " << firstSector + 5 << " dec: " << botLeft << '\n';
0430 
0431     finalTrigger |= (topRight && botLeft);
0432 
0433     ++k;
0434 
0435     if (k > maxTopQuadrants)
0436       break;
0437   }
0438 
0439   return finalTrigger;
0440 }
0441 
0442 // ------------ method called once each job just after ending the event loop  ------------
0443 
0444 void RPCTechnicalTrigger::printinfo() const {
0445   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger::Printing TTU emulators info>" << std::endl;
0446 
0447   for (int k = 0; k < kMaxTtuBoards; ++k) {
0448     m_ttu[k].printinfo();
0449     m_ttuRbcLine[k].printinfo();
0450   }
0451 }
0452 
0453 //define this as a plug-in
0454 DEFINE_FWK_MODULE(RPCTechnicalTrigger);