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