File indexing completed on 2025-06-28 06:19:07
0001 #include "DTHFakeReader.h"
0002 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0003 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0004 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0005 #include "DataFormats/FEDRawData/interface/SLinkRocketHeaders.h"
0006 #include "DataFormats/TCDS/interface/TCDSRaw.h"
0007
0008
0009 #include "EventFilter/Utilities/interface/crc32c.h"
0010
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012
0013 #include "CLHEP/Random/RandGauss.h"
0014
0015 #include <cmath>
0016 #include <sys/time.h>
0017 #include <cstring>
0018 #include <cstdlib>
0019 #include <chrono>
0020
0021
0022 namespace evf {
0023
0024 constexpr unsigned minOrbitBx = 1;
0025 constexpr unsigned maxOrbitBx = 2464;
0026 constexpr unsigned avgEventsPerOrbit = 70;
0027
0028
0029
0030
0031 constexpr unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
0032 constexpr unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);
0033
0034 constexpr double rndFactor = (maxOrbitBx - minOrbitBx + 1) / (double(avgEventsPerOrbit) * RAND_MAX);
0035
0036 DTHFakeReader::DTHFakeReader(const edm::ParameterSet& pset)
0037 : fillRandom_(pset.getUntrackedParameter<bool>("fillRandom", false)),
0038 meansize_(pset.getUntrackedParameter<unsigned int>("meanSize", 1024)),
0039 width_(pset.getUntrackedParameter<unsigned int>("width", 1024)),
0040 injected_errors_per_million_events_(pset.getUntrackedParameter<unsigned int>("injectErrPpm", 0)),
0041 sourceIdList_(
0042 pset.getUntrackedParameter<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>())),
0043 modulo_error_events_(injected_errors_per_million_events_ ? 1000000 / injected_errors_per_million_events_
0044 : 0xffffffff) {
0045 if (fillRandom_) {
0046
0047 auto time_count =
0048 static_cast<long unsigned int>(std::chrono::high_resolution_clock::now().time_since_epoch().count());
0049 std::srand(time_count & 0xffffffff);
0050 }
0051 produces<RawDataBuffer>();
0052 }
0053
0054 void DTHFakeReader::fillRawData(edm::Event& e, RawDataBuffer*& data) {
0055
0056
0057
0058 edm::EventNumber_t orbitId = e.id().event();
0059
0060
0061 std::vector<uint64_t> eventIdList_;
0062 std::map<unsigned int, std::map<uint64_t, uint32_t>> randFedSizes;
0063 for (auto sourceId : sourceIdList_) {
0064 randFedSizes[sourceId] = std::map<uint64_t, uint32_t>();
0065 }
0066 uint32_t totSize = 0;
0067
0068
0069 for (unsigned i = minOrbitBx; i <= maxOrbitBx; i++) {
0070 if ((std::rand() * rndFactor) < 1) {
0071 uint64_t eventId = orbitId * maxOrbitBx + i;
0072 eventIdList_.push_back(eventId);
0073 for (auto sourceId : sourceIdList_) {
0074 float logsiz = CLHEP::RandGauss::shoot(std::log(meansize_), std::log(meansize_) - std::log(width_ / 2.));
0075 size_t size = int(std::exp(logsiz));
0076 size -= size % 16;
0077 if (!size)
0078 size = 16;
0079 randFedSizes[sourceId][eventId] = size;
0080 }
0081 }
0082 }
0083
0084
0085 for (auto sourceId : sourceIdList_) {
0086 auto size = sizeof(DTHOrbitHeader_v1);
0087 for (auto eventId : eventIdList_)
0088 size += randFedSizes[sourceId][eventId] + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
0089 totSize += size;
0090 }
0091 data = new RawDataBuffer(totSize);
0092
0093 for (auto sourceId : sourceIdList_) {
0094 auto size = sizeof(DTHOrbitHeader_v1);
0095 for (auto eventId : eventIdList_)
0096 size += randFedSizes[sourceId][eventId] + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
0097 unsigned char* feddata = data->addSource(sourceId, nullptr, size);
0098
0099 uint64_t fragments_size_bytes = sizeof(DTHOrbitHeader_v1);
0100
0101 uint32_t runningChecksum = 0;
0102 for (auto eventId : eventIdList_) {
0103 unsigned char* fedaddr = feddata + fragments_size_bytes;
0104 fragments_size_bytes +=
0105 fillSLRFED(fedaddr, sourceId, eventId, orbitId, randFedSizes[sourceId][eventId], runningChecksum);
0106 }
0107
0108 new (static_cast<void*>(feddata)) DTHOrbitHeader_v1(sourceId,
0109 e.id().run(),
0110 orbitId,
0111 eventIdList_.size(),
0112 fragments_size_bytes >> evf::DTH_WORD_NUM_BYTES_SHIFT,
0113 0,
0114 runningChecksum);
0115 }
0116 }
0117
0118 void DTHFakeReader::produce(edm::Event& e, edm::EventSetup const& es) {
0119 edm::Handle<RawDataBuffer> rawdata;
0120 RawDataBuffer* fedcoll = nullptr;
0121 fillRawData(e, fedcoll);
0122 std::unique_ptr<RawDataBuffer> bare_product(fedcoll);
0123 e.put(std::move(bare_product));
0124 }
0125
0126 uint32_t DTHFakeReader::fillSLRFED(unsigned char* buf,
0127 const uint32_t sourceId,
0128 edm::EventNumber_t eventId,
0129 const uint32_t orbitId,
0130 uint32_t size,
0131 uint32_t& accum_crc32c) {
0132
0133 const unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
0134 const unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);
0135
0136 uint32_t totsize = size + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
0137 const unsigned fragsize = size + h_size_ + t_size_;
0138
0139
0140 uint8_t emu_status = 2;
0141 uint16_t l1a_types = 1;
0142 uint8_t l1a_phys = 0;
0143 new ((void*)buf) SLinkRocketHeader_v3(sourceId, l1a_types, l1a_phys, emu_status, eventId);
0144
0145
0146 if (fillRandom_) {
0147
0148 size_t size_ui = size - size % sizeof(unsigned int);
0149 for (size_t i = 0; i < size_ui; i += sizeof(unsigned int)) {
0150 *((unsigned int*)(buf + h_size_ + i)) = (unsigned int)std::rand();
0151 }
0152
0153 for (size_t i = size_ui; i < size; i++) {
0154 *(buf + h_size_ + i) = std::rand() & 0xff;
0155 }
0156 }
0157
0158
0159 uint16_t crc = 0;
0160 uint16_t bxid = 0;
0161 uint8_t status = 0;
0162
0163 new ((void*)(buf + h_size_ + size))
0164 SLinkRocketTrailer_v3(status, crc, orbitId, bxid, fragsize >> evf::SLR_WORD_NUM_BYTES_SHIFT, crc);
0165
0166
0167 void* dthTrailerAddr = buf + fragsize;
0168 new (dthTrailerAddr) DTHFragmentTrailer_v1(0, fragsize >> evf::DTH_WORD_NUM_BYTES_SHIFT, eventId, crc);
0169
0170
0171 accum_crc32c = crc32c(accum_crc32c, (const uint8_t*)buf, totsize);
0172
0173 return totsize;
0174 }
0175
0176 void DTHFakeReader::beginLuminosityBlock(edm::LuminosityBlock const& iL, edm::EventSetup const& iE) {
0177 std::cout << "DTHFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
0178 fakeLs_ = iL.luminosityBlock();
0179 }
0180
0181 void DTHFakeReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0182 edm::ParameterSetDescription desc;
0183 desc.setComment("Injector of generated DTH raw orbit fragments for DSAQ testing");
0184 desc.addUntracked<bool>("fillRandom", false);
0185 desc.addUntracked<unsigned int>("meanSize", 1024);
0186 desc.addUntracked<unsigned int>("width", 1024);
0187 desc.addUntracked<unsigned int>("injectErrPpm", 1024);
0188 desc.addUntracked<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>());
0189 descriptions.add("DTHFakeReader", desc);
0190 }
0191 }