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 "HFRingUnpacker.h"
0006 
0007 namespace l1t {
0008   namespace stage1 {
0009     bool HFRingUnpacker::unpack(const Block& block, UnpackerCollections* coll) {
0010       LogDebug("L1T") << "Block ID  = " << block.header().getID() << " size = " << block.header().getSize();
0011 
0012       int nBX = int(ceil(block.header().getSize() / 2.));
0013 
0014       // Find the first and last BXs
0015       int firstBX = -(ceil((double)nBX / 2.) - 1);
0016       int lastBX;
0017       if (nBX % 2 == 0) {
0018         lastBX = ceil((double)nBX / 2.) + 1;
0019       } else {
0020         lastBX = ceil((double)nBX / 2.);
0021       }
0022 
0023       auto resHFBitCounts_ = static_cast<CaloCollections*>(coll)->getCaloSpareHFBitCounts();
0024       resHFBitCounts_->setBXRange(firstBX, lastBX);
0025 
0026       auto resHFRingSums_ = static_cast<CaloCollections*>(coll)->getCaloSpareHFRingSums();
0027       resHFRingSums_->setBXRange(firstBX, lastBX);
0028 
0029       auto reset_ = static_cast<CaloCollections*>(coll)->getEtSums();
0030       reset_->setBXRange(firstBX, lastBX);
0031 
0032       LogDebug("L1T") << "nBX = " << nBX << " first BX = " << firstBX << " lastBX = " << lastBX;
0033 
0034       // Initialise index
0035       int unsigned i = 0;
0036 
0037       // Loop over multiple BX and then number of jets filling jet collection
0038       for (int bx = firstBX; bx < lastBX; bx++) {
0039         uint32_t raw_data0 = block.payload()[i++];
0040         uint32_t raw_data1 = block.payload()[i++];
0041 
0042         /* if (raw_data0 == 0 || raw_data1==0) continue; */
0043 
0044         uint16_t candbit[4];
0045         candbit[0] = raw_data0 & 0xFFFF;
0046         candbit[1] = (raw_data0 >> 16) & 0xFFFF;
0047         candbit[2] = raw_data1 & 0xFFFF;
0048         candbit[3] = (raw_data1 >> 16) & 0xFFFF;
0049 
0050         int hfbitcount = candbit[0] & 0xFFF;
0051         int hfringsum = ((candbit[0] >> 12) & 0x7) | ((candbit[2] & 0x1FF) << 3);
0052         int htmissphi = candbit[1] & 0x1F;
0053         int htmiss = (candbit[1] >> 5) & 0x7F;
0054         int overflowhtmiss = (candbit[1] >> 12) & 0x1;
0055 
0056         l1t::CaloSpare hfbc = l1t::CaloSpare();
0057         hfbc.setHwPt(hfbitcount);
0058         hfbc.setType(l1t::CaloSpare::HFBitCount);
0059         LogDebug("L1T") << "hfbc pT " << hfbc.hwPt();
0060         resHFBitCounts_->push_back(bx, hfbc);
0061 
0062         l1t::CaloSpare hfrs = l1t::CaloSpare();
0063         hfrs.setHwPt(hfringsum);
0064         hfrs.setType(l1t::CaloSpare::HFRingSum);
0065         LogDebug("L1T") << "hfrs pT " << hfrs.hwPt();
0066         resHFRingSums_->push_back(bx, hfrs);
0067 
0068         l1t::EtSum mht{l1t::EtSum::kMissingHt};
0069         mht.setHwPt(htmiss);
0070         mht.setHwPhi(htmissphi);
0071         int flaghtmiss = mht.hwQual();
0072         flaghtmiss |= overflowhtmiss;
0073         mht.setHwQual(flaghtmiss);
0074         LogDebug("L1T") << "MHT: pT " << mht.hwPt() << "is overflow " << overflowhtmiss << std::endl;
0075         reset_->push_back(bx, mht);
0076       }
0077 
0078       return true;
0079     }
0080   }  // namespace stage1
0081 }  // namespace l1t
0082 
0083 DEFINE_L1T_UNPACKER(l1t::stage1::HFRingUnpacker);