File indexing completed on 2024-09-07 04:36:12
0001 #include "EventFilter/L1ScoutingRawToDigi/plugins/ScBMTFRawToDigi.h"
0002
0003 ScBMTFRawToDigi::ScBMTFRawToDigi(const edm::ParameterSet& iConfig) {
0004 using namespace edm;
0005 srcInputTag_ = iConfig.getParameter<InputTag>("srcInputTag");
0006 sourceIdList_ = iConfig.getParameter<std::vector<int>>("sourceIdList");
0007 debug_ = iConfig.getUntrackedParameter<bool>("debug", false);
0008
0009
0010 orbitBuffer_ = std::vector<std::vector<l1ScoutingRun3::BMTFStub>>(3565);
0011 for (auto& bxVec : orbitBuffer_) {
0012 bxVec.reserve(32);
0013 }
0014 nStubsOrbit_ = 0;
0015
0016 produces<l1ScoutingRun3::BMTFStubOrbitCollection>("BMTFStub").setBranchAlias("BMTFStubOrbitCollection");
0017 rawToken_ = consumes<SDSRawDataCollection>(srcInputTag_);
0018 }
0019
0020 ScBMTFRawToDigi::~ScBMTFRawToDigi() {}
0021
0022 void ScBMTFRawToDigi::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0023 using namespace edm;
0024
0025 Handle<SDSRawDataCollection> ScoutingRawDataCollection;
0026 iEvent.getByToken(rawToken_, ScoutingRawDataCollection);
0027
0028 std::unique_ptr<l1ScoutingRun3::BMTFStubOrbitCollection> unpackedStubs(new l1ScoutingRun3::BMTFStubOrbitCollection);
0029
0030 for (const auto& sdsId : sourceIdList_) {
0031 if ((sdsId < SDSNumbering::BmtfMinSDSID) || (sdsId > SDSNumbering::BmtfMaxSDSID))
0032 edm::LogError("ScBMTFRawToDigi::produce")
0033 << "Provided a source ID outside the expected range: " << sdsId << ", expected range ["
0034 << SDSNumbering::BmtfMinSDSID << ", " << SDSNumbering::BmtfMaxSDSID;
0035 const FEDRawData& sourceRawData = ScoutingRawDataCollection->FEDData(sdsId);
0036 size_t orbitSize = sourceRawData.size();
0037
0038 if ((sourceRawData.size() == 0) && debug_) {
0039 std::cout << "No raw data for BMTF FED " << sdsId << std::endl;
0040 }
0041
0042
0043 unpackOrbit(sourceRawData.data(), orbitSize, sdsId);
0044 }
0045
0046
0047 unpackedStubs->fillAndClear(orbitBuffer_, nStubsOrbit_);
0048
0049
0050 iEvent.put(std::move(unpackedStubs), "BMTFStub");
0051 }
0052
0053 void ScBMTFRawToDigi::unpackOrbit(const unsigned char* buf, size_t len, int sdsId) {
0054 using namespace l1ScoutingRun3;
0055
0056
0057 nStubsOrbit_ = 0;
0058
0059 size_t pos = 0;
0060
0061 while (pos < len) {
0062 assert(pos + 4 <= len);
0063
0064 bmtf::block* bl = (bmtf::block*)(buf + pos);
0065
0066 unsigned bx = bl->bx;
0067 unsigned orbit = (bl->orbit) & 0x7FFFFFFF;
0068 unsigned sCount = (bl->header) & 0xff;
0069
0070 size_t pos_increment = 12 + sCount * 8;
0071
0072 assert(pos_increment <= len);
0073
0074 pos += 12;
0075
0076 if (debug_) {
0077 std::cout << " BMTF #" << sdsId << " Orbit " << orbit << ", BX -> " << bx << ", nStubs -> " << sCount
0078 << std::endl;
0079 }
0080
0081
0082 int32_t phi, phiB, tag, qual, eta, qeta, station, wheel, sector;
0083
0084
0085 std::vector<std::vector<bool>> stwh_matrix(4, std::vector<bool>(5, false));
0086 for (unsigned int i = 0; i < sCount; i++) {
0087 uint64_t stub_raw = *(uint64_t*)(buf + pos);
0088 pos += 8;
0089
0090 phi = ((stub_raw >> bmtf::shiftsStubs::phi) & bmtf::masksStubs::phi);
0091 phiB = ((stub_raw >> bmtf::shiftsStubs::phiB) & bmtf::masksStubs::phiB);
0092 qual = ((stub_raw >> bmtf::shiftsStubs::qual) & bmtf::masksStubs::qual);
0093 eta = ((stub_raw >> bmtf::shiftsStubs::eta) & bmtf::masksStubs::eta);
0094 qeta = ((stub_raw >> bmtf::shiftsStubs::qeta) & bmtf::masksStubs::qeta);
0095 station = ((stub_raw >> bmtf::shiftsStubs::station) & bmtf::masksStubs::station) + 1;
0096 wheel = ((stub_raw >> bmtf::shiftsStubs::wheel) & bmtf::masksStubs::wheel);
0097 sector = sdsId - SDSNumbering::BmtfMinSDSID;
0098
0099 if (stwh_matrix[station - 1][wheel + 2] == false) {
0100 tag = 1;
0101 } else {
0102 tag = 0;
0103 }
0104 stwh_matrix[station - 1][wheel + 2] = true;
0105
0106 phi = phi >= 2048 ? phi - 4096 : phi;
0107 phiB = phiB >= 512 ? phiB - 1024 : phiB;
0108 wheel = wheel >= 4 ? wheel - 8 : wheel;
0109
0110 BMTFStub stub(phi, phiB, qual, eta, qeta, station, wheel, sector, tag);
0111 orbitBuffer_[bx].push_back(stub);
0112 nStubsOrbit_++;
0113
0114 if (debug_) {
0115 std::cout << "Stub " << i << ", raw: 0x" << std::hex << stub_raw << std::dec << std::endl;
0116 std::cout << "\tPhi: " << phi << std::endl;
0117 std::cout << "\tPhiB: " << phiB << std::endl;
0118 std::cout << "\tQuality: " << qual << std::endl;
0119 std::cout << "\tEta: " << eta << std::endl;
0120 std::cout << "\tQEta: " << qeta << std::endl;
0121 std::cout << "\tStation: " << station << std::endl;
0122 std::cout << "\tWheel: " << wheel << std::endl;
0123 std::cout << "\tSector: " << sector << std::endl;
0124 std::cout << "\tTag: " << tag << std::endl;
0125 }
0126 }
0127
0128 }
0129 }
0130
0131 void ScBMTFRawToDigi::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0132 edm::ParameterSetDescription desc;
0133 desc.setUnknown();
0134 descriptions.addDefault(desc);
0135 }
0136
0137 DEFINE_FWK_MODULE(ScBMTFRawToDigi);