File indexing completed on 2024-04-06 12:10:51
0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "EventFilter/L1TRawToDigi/plugins/UnpackerFactory.h"
0003
0004 #include "GTCollections.h"
0005 #include "GlobalAlgBlkUnpacker.h"
0006
0007 namespace l1t {
0008 namespace stage2 {
0009 bool GlobalAlgBlkUnpacker::unpack(const Block& block, UnpackerCollections* coll) {
0010 LogDebug("L1T") << "AMCNo " << block.amc().getAMCNumber() << " Block ID = " << block.header().getID()
0011 << " size = " << block.header().getSize();
0012
0013
0014
0015 unsigned int wdPerBX = 6;
0016 unsigned int initialBlkID = 33;
0017 unsigned int intermBlkID = 39;
0018 unsigned int finalBlkID = 45;
0019
0020
0021 unsigned int uGTBoard = block.amc().getAMCNumber() - 1;
0022
0023 int nBX =
0024 int(ceil(block.header().getSize() / 6.));
0025
0026
0027 int firstBX = -(ceil((double)nBX / 2.) - 1);
0028 int lastBX;
0029 if (nBX % 2 == 0) {
0030 lastBX = ceil((double)nBX / 2.);
0031 } else {
0032 lastBX = ceil((double)nBX / 2.) - 1;
0033 }
0034
0035 auto res_ = static_cast<GTCollections*>(coll)->getAlgs();
0036 res_->setBXRange(firstBX, lastBX);
0037
0038 LogDebug("L1T") << "nBX = " << nBX << " first BX = " << firstBX << " lastBX = " << lastBX << endl;
0039
0040
0041 int numBX = 0;
0042 for (int bx = firstBX; bx <= lastBX; bx++) {
0043
0044 if (block.header().getID() == initialBlkID && uGTBoard == 0) {
0045 LogDebug("L1T") << "Creating GT Algorithm Block for BX =" << bx << std::endl;
0046 GlobalAlgBlk talg = GlobalAlgBlk();
0047 res_->push_back(bx, talg);
0048 }
0049
0050 else if (res_->isEmpty(bx))
0051 throw cms::Exception("InvalidGlobalAlgBlkBxCollection")
0052 << "The GlobalAlgBlk unpacker result vector is empty, but is not receiving the first expected header "
0053 "ID! This may be due to corrupted, or poorly formatted events.\n"
0054 << "uGTBoard: " << uGTBoard << "\nBX: " << bx << "\nFirst expected block: " << initialBlkID
0055 << "\nReceived block: " << block.header().getID();
0056
0057
0058 GlobalAlgBlk alg = res_->at(bx, 0);
0059
0060
0061
0062
0063
0064 int algOffset = (block.header().getID() - initialBlkID + 1) / 2;
0065 algOffset = (algOffset % 3) * 192;
0066
0067 for (unsigned int wd = 0; wd < wdPerBX; wd++) {
0068 uint32_t raw_data = block.payload()[wd + numBX * wdPerBX];
0069 LogDebug("L1T") << "BX " << bx << " payload word " << wd << " 0x" << hex << raw_data << " offset=" << dec
0070 << algOffset << std::endl;
0071
0072
0073 if ((block.header().getID() != initialBlkID + 4 && block.header().getID() != intermBlkID + 4 &&
0074 block.header().getID() != finalBlkID + 4) ||
0075 wd < 4) {
0076 for (unsigned int bt = 0; bt < 32; bt++) {
0077 int val = ((raw_data >> bt) & 0x1);
0078 unsigned int algBit = bt + wd * 32 + algOffset;
0079
0080 if (val == 1 && algBit < alg.maxPhysicsTriggers) {
0081 LogDebug("L1T") << "Found valid alg bit (" << algBit << ") on bit (" << bt << ") word (" << wd
0082 << ") algOffset (" << algOffset << ") block ID (" << block.header().getID() << ")"
0083 << " Board# " << uGTBoard << std::endl;
0084 if (block.header().getID() < initialBlkID + 5) {
0085 alg.setAlgoDecisionInitial(algBit, true);
0086 } else if (block.header().getID() < intermBlkID + 5) {
0087 alg.setAlgoDecisionInterm(algBit, true);
0088 } else {
0089 alg.setAlgoDecisionFinal(algBit, true);
0090 }
0091 } else if (val == 1) {
0092 LogDebug("L1T") << "Found invalid alg bit (" << algBit << ") out of range on bit (" << bt << ") word ("
0093 << wd << ") algOffset (" << algOffset << ") block ID (" << block.header().getID() << ")"
0094 << std::endl;
0095 }
0096 }
0097
0098 } else if (block.header().getID() == initialBlkID + 4 && (wd == 4 || wd == 5)) {
0099
0100 if (wd == 4)
0101 alg.setL1MenuUUID(raw_data);
0102
0103 if (wd == 5)
0104 alg.setL1FirmwareUUID(raw_data);
0105
0106 } else if (block.header().getID() == finalBlkID + 4 && wd == 4) {
0107
0108 if ((raw_data & 0x100) >> 8)
0109 alg.setFinalORVeto(true);
0110 if ((raw_data & 0x1) >> 0)
0111 alg.setFinalORPreVeto(true);
0112 LogDebug("L1T") << " Packing the FinalOR " << wd << " 0x" << hex << raw_data << endl;
0113 } else if (block.header().getID() == finalBlkID + 4 && wd == 5) {
0114
0115 alg.setPreScColumn(raw_data & 0xFF);
0116 LogDebug("L1T") << " Packing the Prescale Column " << wd << " 0x" << hex << raw_data << endl;
0117 }
0118 }
0119
0120
0121
0122 alg.setFinalOR(alg.getFinalORPreVeto() && !alg.getFinalORVeto());
0123
0124
0125 res_->set(bx, 0, alg);
0126
0127
0128
0129
0130 numBX++;
0131 }
0132
0133 return true;
0134 }
0135 }
0136 }
0137
0138 DEFINE_L1T_UNPACKER(l1t::stage2::GlobalAlgBlkUnpacker);