Back to home page

Project CMSSW displayed by LXR

 
 

    


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 "L1Trigger/L1TCalorimeter/interface/CaloTools.h"
0005 
0006 #include "L1TObjectCollections.h"
0007 
0008 #include "L1TStage2Layer2Constants.h"
0009 #include "JetUnpacker.h"
0010 
0011 namespace l1t {
0012   namespace stage2 {
0013     JetUnpacker::JetUnpacker() : JetCopy_(0) {}
0014 
0015     bool JetUnpacker::unpack(const Block& block, UnpackerCollections* coll) {
0016       using namespace l1t::stage2::layer2;
0017 
0018       LogDebug("L1T") << "Block ID  = " << block.header().getID() << " size = " << block.header().getSize();
0019 
0020       int nBX = int(ceil(block.header().getSize() / (double)demux::nOutputFramePerBX));  // 6 frames per BX
0021 
0022       // Find the first and last BXs
0023       int firstBX = -(ceil((double)nBX / 2.) - 1);
0024       int lastBX;
0025       if (nBX % 2 == 0) {
0026         lastBX = ceil((double)nBX / 2.);
0027       } else {
0028         lastBX = ceil((double)nBX / 2.) - 1;
0029       }
0030 
0031       auto res_ = static_cast<L1TObjectCollections*>(coll)->getJets(JetCopy_);
0032       res_->setBXRange(firstBX, lastBX);
0033 
0034       LogDebug("L1T") << "nBX = " << nBX << " first BX = " << firstBX << " lastBX = " << lastBX;
0035 
0036       // Loop over multiple BX and then number of jets filling jet collection
0037       for (int bx = firstBX; bx <= lastBX; bx++) {
0038         for (unsigned iJet = 0; iJet < demux::nJetPerLink && iJet < block.header().getSize(); iJet++) {
0039           int iFrame = (bx - firstBX) * demux::nOutputFramePerBX + iJet;
0040           uint32_t raw_data = block.payload().at(iFrame);
0041 
0042           if (raw_data == 0)
0043             continue;
0044 
0045           l1t::Jet jet = l1t::Jet();
0046 
0047           jet.setHwPt(raw_data & 0x7FF);
0048 
0049           if (jet.hwPt() == 0)
0050             continue;
0051 
0052           int abs_eta = (raw_data >> 11) & 0x7F;
0053           if ((raw_data >> 18) & 0x1) {
0054             jet.setHwEta(-1 * (128 - abs_eta));
0055           } else {
0056             jet.setHwEta(abs_eta);
0057           }
0058 
0059           jet.setHwPhi((raw_data >> 19) & 0xFF);
0060           jet.setHwQual((raw_data >> 27) & 0x7);  // Assume 3 bits for now? Leaves 2 bits spare
0061 
0062           LogDebug("L1T") << "Jet: eta " << jet.hwEta() << " phi " << jet.hwPhi() << " pT " << jet.hwPt() << " qual "
0063                           << jet.hwQual() << " bx " << bx;
0064 
0065           jet.setP4(l1t::CaloTools::p4Demux(&jet));
0066 
0067           res_->push_back(bx, jet);
0068         }
0069       }
0070 
0071       return true;
0072     }
0073   }  // namespace stage2
0074 }  // namespace l1t
0075 
0076 DEFINE_L1T_UNPACKER(l1t::stage2::JetUnpacker);