File indexing completed on 2024-04-06 12:11:06
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
0046 if (tcdsFEDID_ < FEDNumbering::MINTCDSuTCAFEDID)
0047 throw cms::Exception("DaqFakeReader::DaqFakeReader")
0048 << " TCDS FED ID lower than " << FEDNumbering::MINTCDSuTCAFEDID;
0049 if (fillRandom_) {
0050
0051 auto time_count =
0052 static_cast<long unsigned int>(std::chrono::high_resolution_clock::now().time_since_epoch().count());
0053 srand(time_count & 0xffffffff);
0054 }
0055 produces<FEDRawDataCollection>();
0056 }
0057
0058
0059 DaqFakeReader::~DaqFakeReader() {}
0060
0061
0062
0063
0064
0065
0066 int DaqFakeReader::fillRawData(Event& e, FEDRawDataCollection*& data) {
0067
0068 data = new FEDRawDataCollection();
0069 EventID eID = e.id();
0070 auto ls = e.luminosityBlock();
0071
0072 if (!empty_events) {
0073
0074 eventNum++;
0075
0076
0077 fillFEDs(FEDNumbering::MINSiPixelFEDID, FEDNumbering::MAXSiPixelFEDID, eID, *data, meansize, width);
0078 fillFEDs(FEDNumbering::MINSiStripFEDID, FEDNumbering::MAXSiStripFEDID, eID, *data, meansize, width);
0079 fillFEDs(FEDNumbering::MINDTFEDID, FEDNumbering::MAXDTFEDID, eID, *data, meansize, width);
0080 fillFEDs(FEDNumbering::MINCSCFEDID, FEDNumbering::MAXCSCFEDID, eID, *data, meansize, width);
0081 fillFEDs(FEDNumbering::MINRPCFEDID, FEDNumbering::MAXRPCFEDID, eID, *data, meansize, width);
0082 fillFEDs(FEDNumbering::MINECALFEDID, FEDNumbering::MAXECALFEDID, eID, *data, meansize, width);
0083 fillFEDs(FEDNumbering::MINHCALFEDID, FEDNumbering::MAXHCALFEDID, eID, *data, meansize, width);
0084
0085 timeval now;
0086 gettimeofday(&now, nullptr);
0087 fillTCDSFED(eID, *data, ls, &now);
0088 }
0089 return 1;
0090 }
0091
0092 void DaqFakeReader::produce(Event& e, EventSetup const& es) {
0093 edm::Handle<FEDRawDataCollection> rawdata;
0094 FEDRawDataCollection* fedcoll = nullptr;
0095 fillRawData(e, fedcoll);
0096 std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0097 e.put(std::move(bare_product));
0098 }
0099
0100
0101 void DaqFakeReader::fillFEDs(
0102 const int fedmin, const int fedmax, EventID& eID, FEDRawDataCollection& data, float meansize, float width) {
0103
0104 for (int fedId = fedmin; fedId <= fedmax; ++fedId) {
0105
0106 float logsiz = CLHEP::RandGauss::shoot(std::log(meansize), std::log(meansize) - std::log(width / 2.));
0107 size_t size = int(std::exp(logsiz));
0108 size -= size % 8;
0109
0110 FEDRawData& feddata = data.FEDData(fedId);
0111
0112 feddata.resize(size + 16);
0113
0114 if (fillRandom_) {
0115
0116 size_t size_ui = size - size % sizeof(unsigned int);
0117 for (size_t i = 0; i < size_ui; i += sizeof(unsigned int)) {
0118 *((unsigned int*)(feddata.data() + i)) = (unsigned int)rand();
0119 }
0120
0121 for (size_t i = size_ui; i < size; i++) {
0122 *(feddata.data() + i) = rand() & 0xff;
0123 }
0124 }
0125
0126
0127 FEDHeader::set(feddata.data(),
0128 1,
0129 eID.event(),
0130 0,
0131 fedId);
0132
0133
0134
0135
0136 int crc = 0;
0137 FEDTrailer::set(feddata.data() + 8 + size,
0138 size / 8 + 2,
0139 crc,
0140 0,
0141 0);
0142 }
0143 }
0144
0145 void DaqFakeReader::fillTCDSFED(EventID& eID, FEDRawDataCollection& data, uint32_t ls, timeval* now) {
0146 uint32_t fedId = tcdsFEDID_;
0147 FEDRawData& feddata = data.FEDData(fedId);
0148 uint32_t size = sizeof(tcds::Raw_v1);
0149 feddata.resize(size + 16);
0150
0151 uint64_t orbitnr = 0;
0152 uint16_t bxid = 0;
0153
0154 FEDHeader::set(feddata.data(),
0155 1,
0156 eID.event(),
0157 bxid,
0158 fedId);
0159
0160 tcds::Raw_v1* tcds = reinterpret_cast<tcds::Raw_v1*>(feddata.data() + FEDHeader::length);
0161 tcds::BST_v1* bst = const_cast<tcds::BST_v1*>(&tcds->bst);
0162 tcds::Header_v1* header = const_cast<tcds::Header_v1*>(&tcds->header);
0163
0164 const_cast<uint32_t&>(bst->gpstimehigh) = now->tv_sec;
0165 const_cast<uint32_t&>(bst->gpstimelow) = now->tv_usec;
0166 const_cast<uint16_t&>(bst->lhcFillHigh) = 0;
0167 const_cast<uint16_t&>(bst->lhcFillLow) = 0;
0168
0169 const_cast<uint32_t&>(header->orbitHigh) = orbitnr & 0xffff00;
0170 const_cast<uint16_t&>(header->orbitLow) = orbitnr & 0xff;
0171 const_cast<uint16_t&>(header->bxid) = bxid;
0172
0173 const_cast<uint64_t&>(header->eventNumber) = eID.event();
0174 const_cast<uint32_t&>(header->lumiSection) = ls;
0175
0176 int crc = 0;
0177 FEDTrailer::set(feddata.data() + 8 + size,
0178 size / 8 + 2,
0179 crc,
0180 0,
0181 0);
0182 }
0183
0184 void DaqFakeReader::beginLuminosityBlock(LuminosityBlock const& iL, EventSetup const& iE) {
0185 std::cout << "DaqFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
0186 fakeLs_ = iL.luminosityBlock();
0187 }
0188
0189 void DaqFakeReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0190 edm::ParameterSetDescription desc;
0191 desc.setComment("Injector of generated raw FED data for DAQ testing");
0192 desc.addUntracked<bool>("emptyEvents", false);
0193 desc.addUntracked<bool>("fillRandom", false);
0194 desc.addUntracked<unsigned int>("meanSize", 1024);
0195 desc.addUntracked<unsigned int>("width", 1024);
0196 desc.addUntracked<unsigned int>("injectErrPpm", 1024);
0197 desc.addUntracked<unsigned int>("tcdsFEDID", 1024);
0198 descriptions.add("DaqFakeReader", desc);
0199 }