File indexing completed on 2025-06-29 22:57:58
0001
0002
0003
0004
0005
0006 #include "DaqFakeReader.h"
0007 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0008 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0009 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0010 #include "DataFormats/TCDS/interface/TCDSRaw.h"
0011
0012 #include "EventFilter/Utilities/interface/GlobalEventNumber.h"
0013
0014 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0015
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017
0018 #include "CLHEP/Random/RandGauss.h"
0019
0020 #include <cmath>
0021 #include <sys/time.h>
0022 #include <cstring>
0023 #include <cstdlib>
0024 #include <chrono>
0025
0026 using namespace std;
0027 using namespace edm;
0028
0029
0030
0031
0032
0033
0034 DaqFakeReader::DaqFakeReader(const edm::ParameterSet& pset)
0035 : runNum(1),
0036 eventNum(1),
0037 empty_events(pset.getUntrackedParameter<bool>("emptyEvents", false)),
0038 fillRandom_(pset.getUntrackedParameter<bool>("fillRandom", false)),
0039 meansize(pset.getUntrackedParameter<unsigned int>("meanSize", 1024)),
0040 width(pset.getUntrackedParameter<unsigned int>("width", 1024)),
0041 injected_errors_per_million_events(pset.getUntrackedParameter<unsigned int>("injectErrPpm", 0)),
0042 tcdsFEDID_(pset.getUntrackedParameter<unsigned int>("tcdsFEDID", 1024)),
0043 modulo_error_events(injected_errors_per_million_events ? 1000000 / injected_errors_per_million_events
0044 : 0xffffffff),
0045 subsystems_(pset.getUntrackedParameter<std::vector<std::string>>("subsystems")) {
0046 for (auto const& subsystem : subsystems_) {
0047 if (subsystem == "TCDS")
0048 haveTCDS_ = true;
0049 else if (subsystem == "SiPixel")
0050 haveSiPixel_ = true;
0051 else if (subsystem == "SiStrip")
0052 haveSiStrip_ = true;
0053 else if (subsystem == "ECAL")
0054 haveECAL_ = true;
0055 else if (subsystem == "HCAL")
0056 haveHCAL_ = true;
0057 else if (subsystem == "DT")
0058 haveDT_ = true;
0059 else if (subsystem == "CSC")
0060 haveCSC_ = true;
0061 else if (subsystem == "RPC")
0062 haveRPC_ = true;
0063 }
0064
0065 if (haveTCDS_ && tcdsFEDID_ < FEDNumbering::MINTCDSuTCAFEDID)
0066 throw cms::Exception("DaqFakeReader::DaqFakeReader")
0067 << " TCDS FED ID lower than " << FEDNumbering::MINTCDSuTCAFEDID;
0068 if (fillRandom_) {
0069
0070 auto time_count =
0071 static_cast<long unsigned int>(std::chrono::high_resolution_clock::now().time_since_epoch().count());
0072 srand(time_count & 0xffffffff);
0073 }
0074 produces<FEDRawDataCollection>();
0075 }
0076
0077
0078 DaqFakeReader::~DaqFakeReader() {}
0079
0080
0081
0082
0083
0084
0085 int DaqFakeReader::fillRawData(Event& e, FEDRawDataCollection*& data) {
0086
0087 data = new FEDRawDataCollection();
0088 EventID eID = e.id();
0089 auto ls = e.luminosityBlock();
0090
0091 if (!empty_events) {
0092
0093 eventNum++;
0094
0095
0096 if (haveSiPixel_)
0097 fillFEDs(FEDNumbering::MINSiPixelFEDID, FEDNumbering::MAXSiPixelFEDID, eID, *data, meansize, width);
0098 if (haveSiStrip_)
0099 fillFEDs(FEDNumbering::MINSiStripFEDID, FEDNumbering::MAXSiStripFEDID, eID, *data, meansize, width);
0100 if (haveECAL_)
0101 fillFEDs(FEDNumbering::MINECALFEDID, FEDNumbering::MAXECALFEDID, eID, *data, meansize, width);
0102 if (haveHCAL_)
0103 fillFEDs(FEDNumbering::MINHCALFEDID, FEDNumbering::MAXHCALFEDID, eID, *data, meansize, width);
0104 if (haveDT_)
0105 fillFEDs(FEDNumbering::MINDTFEDID, FEDNumbering::MAXDTFEDID, eID, *data, meansize, width);
0106 if (haveCSC_)
0107 fillFEDs(FEDNumbering::MINCSCFEDID, FEDNumbering::MAXCSCFEDID, eID, *data, meansize, width);
0108 if (haveRPC_)
0109 fillFEDs(FEDNumbering::MINRPCFEDID, FEDNumbering::MAXRPCFEDID, eID, *data, meansize, width);
0110
0111 timeval now;
0112 gettimeofday(&now, nullptr);
0113 if (haveTCDS_)
0114 fillTCDSFED(eID, *data, ls, &now);
0115 }
0116 return 1;
0117 }
0118
0119 void DaqFakeReader::produce(Event& e, EventSetup const& es) {
0120 edm::Handle<FEDRawDataCollection> rawdata;
0121 FEDRawDataCollection* fedcoll = nullptr;
0122 fillRawData(e, fedcoll);
0123 std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0124 e.put(std::move(bare_product));
0125 }
0126
0127
0128 void DaqFakeReader::fillFEDs(
0129 const int fedmin, const int fedmax, EventID& eID, FEDRawDataCollection& data, float meansize, float width) {
0130
0131 for (int fedId = fedmin; fedId <= fedmax; ++fedId) {
0132
0133 float logsiz = CLHEP::RandGauss::shoot(std::log(meansize), std::log(meansize) - std::log(width / 2.));
0134 size_t size = int(std::exp(logsiz));
0135 size -= size % 8;
0136
0137 FEDRawData& feddata = data.FEDData(fedId);
0138
0139 feddata.resize(size + 16);
0140
0141 if (fillRandom_) {
0142
0143 size_t size_ui = size - size % sizeof(unsigned int);
0144 for (size_t i = 0; i < size_ui; i += sizeof(unsigned int)) {
0145 *((unsigned int*)(feddata.data() + i)) = (unsigned int)rand();
0146 }
0147
0148 for (size_t i = size_ui; i < size; i++) {
0149 *(feddata.data() + i) = rand() & 0xff;
0150 }
0151 }
0152
0153
0154 FEDHeader::set(feddata.data(),
0155 1,
0156 eID.event(),
0157 0,
0158 fedId);
0159
0160
0161
0162
0163 int crc = 0;
0164 FEDTrailer::set(feddata.data() + 8 + size,
0165 size / 8 + 2,
0166 crc,
0167 0,
0168 0);
0169 }
0170 }
0171
0172 void DaqFakeReader::fillTCDSFED(EventID& eID, FEDRawDataCollection& data, uint32_t ls, timeval* now) {
0173 uint32_t fedId = tcdsFEDID_;
0174 FEDRawData& feddata = data.FEDData(fedId);
0175 uint32_t size = sizeof(tcds::Raw_v1);
0176 feddata.resize(size + 16);
0177
0178 uint64_t orbitnr = 0;
0179 uint16_t bxid = 0;
0180
0181 FEDHeader::set(feddata.data(),
0182 1,
0183 eID.event(),
0184 bxid,
0185 fedId);
0186
0187 tcds::Raw_v1* tcds = reinterpret_cast<tcds::Raw_v1*>(feddata.data() + FEDHeader::length);
0188 tcds::BST_v1* bst = const_cast<tcds::BST_v1*>(&tcds->bst);
0189 tcds::Header_v1* header = const_cast<tcds::Header_v1*>(&tcds->header);
0190
0191 const_cast<uint32_t&>(bst->gpstimehigh) = now->tv_sec;
0192 const_cast<uint32_t&>(bst->gpstimelow) = now->tv_usec;
0193 const_cast<uint16_t&>(bst->lhcFillHigh) = 0;
0194 const_cast<uint16_t&>(bst->lhcFillLow) = 0;
0195
0196 const_cast<uint32_t&>(header->orbitHigh) = orbitnr & 0xffff00;
0197 const_cast<uint16_t&>(header->orbitLow) = orbitnr & 0xff;
0198 const_cast<uint16_t&>(header->bxid) = bxid;
0199
0200 const_cast<uint64_t&>(header->eventNumber) = eID.event();
0201 const_cast<uint32_t&>(header->lumiSection) = ls;
0202
0203 int crc = 0;
0204 FEDTrailer::set(feddata.data() + 8 + size,
0205 size / 8 + 2,
0206 crc,
0207 0,
0208 0);
0209 }
0210
0211 void DaqFakeReader::beginLuminosityBlock(LuminosityBlock const& iL, EventSetup const& iE) {
0212 std::cout << "DaqFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
0213 fakeLs_ = iL.luminosityBlock();
0214 }
0215
0216 void DaqFakeReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0217 edm::ParameterSetDescription desc;
0218 desc.setComment("Injector of generated raw FED data for DAQ testing");
0219 desc.addUntracked<bool>("emptyEvents", false);
0220 desc.addUntracked<bool>("fillRandom", false);
0221 desc.addUntracked<unsigned int>("meanSize", 1024);
0222 desc.addUntracked<unsigned int>("width", 1024);
0223 desc.addUntracked<unsigned int>("injectErrPpm", 1024);
0224 desc.addUntracked<unsigned int>("tcdsFEDID", 1024);
0225 desc.addUntracked<std::vector<std::string>>(
0226 "subsystems",
0227 std::initializer_list<std::string>({"TCDS", "SiPixel", "SiStrip", "ECAL", "HCAL", "DT", "CSC", "RPC"}));
0228 descriptions.add("DaqFakeReader", desc);
0229 }