File indexing completed on 2024-04-06 12:10:50
0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "EventFilter/L1TRawToDigi/plugins/UnpackerFactory.h"
0003
0004 #include "L1Trigger/L1TCalorimeter/interface/CaloTools.h"
0005
0006 #include "CaloCollections.h"
0007 #include "CaloTowerUnpacker.h"
0008 #include "L1TStage2Layer2Constants.h"
0009
0010 namespace l1t {
0011 namespace stage2 {
0012 bool CaloTowerUnpacker::unpack(const Block& block, UnpackerCollections* coll) {
0013
0014 unsigned int tmt = block.amc().getBoardID() - l1t::stage2::layer2::mp::offsetBoardId + 1;
0015 unsigned int bxid = block.amc().getBX();
0016
0017 LogDebug("L1T") << "Unpacking TMT # " << tmt << " for BX " << bxid;
0018
0019
0020 unsigned link = block.header().getID() / 2;
0021
0022
0023 unsigned link_phi = (link % 2 == 0) ? link : (link - 1);
0024
0025
0026 unsigned nframes = 40;
0027
0028 int nBX = int(ceil(
0029 block.header().getSize() /
0030 nframes));
0031
0032
0033 int firstBX = -(std::ceil((double)nBX / 2.) - 1);
0034 int lastBX;
0035 if (nBX % 2 == 0) {
0036 lastBX = std::ceil((double)nBX / 2.);
0037 } else {
0038 lastBX = std::ceil((double)nBX / 2.) - 1;
0039 }
0040
0041 auto res_ = static_cast<CaloCollections*>(coll)->getTowers();
0042 res_->setBXRange(std::min(firstBX, res_->getFirstBX()), std::max(lastBX, res_->getLastBX()));
0043
0044 LogDebug("L1T") << "Block : id=" << block.header().getID() << ", size=" << block.header().getSize()
0045 << ", link=" << link << ", link_phi=" << link_phi << ", nBX=" << nBX << ", firstBX=" << firstBX
0046 << ", lastBX=" << lastBX;
0047
0048
0049 for (int bx = firstBX; bx <= lastBX; bx++) {
0050 for (unsigned iframe = 0; iframe < nframes && iframe < block.header().getSize(); ++iframe) {
0051 uint32_t raw_data = block.payload().at(iframe);
0052
0053 if ((raw_data & 0xFFFF) != 0) {
0054 l1t::CaloTower tower1 = l1t::CaloTower();
0055
0056
0057 tower1.setHwPt(raw_data & 0x1FF);
0058 tower1.setHwQual((raw_data >> 12) & 0xF);
0059 tower1.setHwEtRatio((raw_data >> 9) & 0x7);
0060 tower1.setHwPhi(link_phi + 1);
0061
0062 int ieta = iframe + 1;
0063 if (link % 2 != 0)
0064 ieta = ieta * -1;
0065
0066 tower1.setHwEta(CaloTools::caloEta(ieta));
0067
0068 LogDebug("L1T") << "Tower 1: Eta " << tower1.hwEta() << " phi " << tower1.hwPhi() << " pT " << tower1.hwPt()
0069 << " frame " << iframe << " qual " << tower1.hwQual() << " EtRatio " << tower1.hwEtRatio();
0070
0071 res_->push_back(bx, tower1);
0072 }
0073
0074 if (((raw_data >> 16) & 0xFFFF) != 0) {
0075
0076 l1t::CaloTower tower2 = l1t::CaloTower();
0077
0078 tower2.setHwPt((raw_data >> 16) & 0x1FF);
0079 tower2.setHwQual((raw_data >> 28) & 0xF);
0080 tower2.setHwEtRatio((raw_data >> 25) & 0x7);
0081 tower2.setHwPhi(link_phi + 2);
0082
0083 int ieta = iframe + 1;
0084 if (link % 2 != 0)
0085 ieta = ieta * -1;
0086 tower2.setHwEta(CaloTools::caloEta(ieta));
0087
0088 LogDebug("L1T") << "Tower 2: Eta " << tower2.hwEta() << " phi " << tower2.hwPhi() << " pT " << tower2.hwPt()
0089 << " frame " << iframe << " qual " << tower2.hwQual() << " EtRatio " << tower2.hwEtRatio();
0090
0091 res_->push_back(bx, tower2);
0092 }
0093 }
0094 }
0095
0096 return true;
0097 }
0098 }
0099 }
0100
0101 DEFINE_L1T_UNPACKER(l1t::stage2::CaloTowerUnpacker);