Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:48

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "EventFilter/L1TRawToDigi/plugins/UnpackerFactory.h"
0003 
0004 #include "CaloCollections.h"
0005 #include "CaloSpareHFUnpacker.h"
0006 
0007 namespace l1t {
0008   namespace stage1 {
0009     bool CaloSpareHFUnpacker::unpack(const Block& block, UnpackerCollections* coll) {
0010       LogDebug("L1T") << "Block ID  = " << block.header().getID() << " size = " << block.header().getSize();
0011 
0012       int nBX, firstBX, lastBX;
0013       nBX = int(ceil(block.header().getSize() / 2.));
0014       getBXRange(nBX, firstBX, lastBX);
0015 
0016       auto resHFBitCounts_ = static_cast<CaloCollections*>(coll)->getCaloSpareHFBitCounts();
0017       resHFBitCounts_->setBXRange(firstBX, lastBX);
0018 
0019       auto resHFRingSums_ = static_cast<CaloCollections*>(coll)->getCaloSpareHFRingSums();
0020       resHFRingSums_->setBXRange(firstBX, lastBX);
0021 
0022       LogDebug("L1T") << "nBX = " << nBX << " first BX = " << firstBX << " lastBX = " << lastBX;
0023 
0024       // Initialise index
0025       int unsigned i = 0;
0026 
0027       // Loop over multiple BX and then number of jets filling jet collection
0028       for (int bx = firstBX; bx <= lastBX; bx++) {
0029         uint32_t raw_data0 = block.payload()[i++];
0030         uint32_t raw_data1 = block.payload()[i++];
0031 
0032         /* if (raw_data0 == 0 || raw_data1==0) continue; */
0033 
0034         uint16_t candbit[2];
0035         candbit[0] = raw_data0 & 0xFFFF;
0036         candbit[1] = raw_data1 & 0xFFFF;
0037 
0038         int hfbitcount = candbit[0] & 0xFFF;
0039         int hfringsum = ((candbit[0] >> 12) & 0x7) | ((candbit[1] & 0x1FF) << 3);
0040 
0041         l1t::CaloSpare hfbc = l1t::CaloSpare();
0042         hfbc.setHwPt(hfbitcount);
0043         hfbc.setType(l1t::CaloSpare::HFBitCount);
0044         LogDebug("L1T") << "hfbc pT " << hfbc.hwPt();
0045         resHFBitCounts_->push_back(bx, hfbc);
0046 
0047         l1t::CaloSpare hfrs = l1t::CaloSpare();
0048         hfrs.setHwPt(hfringsum);
0049         hfrs.setType(l1t::CaloSpare::HFRingSum);
0050         LogDebug("L1T") << "hfrs pT " << hfrs.hwPt();
0051         resHFRingSums_->push_back(bx, hfrs);
0052       }
0053 
0054       return true;
0055     }
0056   }  // namespace stage1
0057 }  // namespace l1t
0058 
0059 DEFINE_L1T_UNPACKER(l1t::stage1::CaloSpareHFUnpacker);