Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-12 23:41:43

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/TCDS/interface/TCDSRaw.h"
0006 
0007 //#include "EventFilter/Utilities/interface/GlobalEventNumber.h"
0008 #include "EventFilter/Utilities/interface/crc32c.h"
0009 
0010 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0011 
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 
0014 #include "CLHEP/Random/RandGauss.h"
0015 
0016 #include <cmath>
0017 #include <sys/time.h>
0018 #include <cstring>
0019 #include <cstdlib>
0020 #include <chrono>
0021 
0022 //using namespace edm;
0023 namespace evf {
0024 
0025   constexpr unsigned minOrbitBx = 1;
0026   constexpr unsigned maxOrbitBx = 2464;
0027   constexpr unsigned avgEventsPerOrbit = 70;
0028 
0029   //constexpr unsigned h_size_ = 8;//for SLink FEDs
0030   //constexpr unsigned t_size_ = 8;
0031 
0032   constexpr unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
0033   constexpr unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);
0034 
0035   constexpr double rndFactor = (maxOrbitBx - minOrbitBx + 1) / (double(avgEventsPerOrbit) * RAND_MAX);
0036 
0037   DTHFakeReader::DTHFakeReader(const edm::ParameterSet& pset)
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         sourceIdList_(
0043             pset.getUntrackedParameter<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>())),
0044         modulo_error_events_(injected_errors_per_million_events_ ? 1000000 / injected_errors_per_million_events_
0045                                                                  : 0xffffffff) {
0046     if (fillRandom_) {
0047       //intialize random seed
0048       auto time_count =
0049           static_cast<long unsigned int>(std::chrono::high_resolution_clock::now().time_since_epoch().count());
0050       std::srand(time_count & 0xffffffff);
0051     }
0052     produces<FEDRawDataCollection>();
0053   }
0054 
0055   void DTHFakeReader::fillRawData(edm::Event& e, FEDRawDataCollection*& data) {
0056     // a null pointer is passed, need to allocate the fed collection (reusing it as container)
0057     data = new FEDRawDataCollection();
0058     //auto ls = e.luminosityBlock();
0059     //this will be used as orbit counter
0060     edm::EventNumber_t orbitId = e.id().event();
0061 
0062     //generate eventID. Orbits start from 0 or 1?
0063     std::vector<uint64_t> eventIdList_;
0064     std::map<unsigned int, std::map<uint64_t, uint32_t>> randFedSizes;
0065     for (auto sourceId : sourceIdList_) {
0066       randFedSizes[sourceId] = std::map<uint64_t, uint32_t>();
0067     }
0068 
0069     //randomize which orbit was accepted
0070     for (unsigned i = minOrbitBx; i <= maxOrbitBx; i++) {
0071       if ((std::rand() * rndFactor) < 1) {
0072         uint64_t eventId = orbitId * maxOrbitBx + i;
0073         eventIdList_.push_back(eventId);
0074         for (auto sourceId : sourceIdList_) {
0075           float logsiz = CLHEP::RandGauss::shoot(std::log(meansize_), std::log(meansize_) - std::log(width_ / 2.));
0076           size_t size = int(std::exp(logsiz));
0077           size -= size % 16;  // all blocks aligned to 128 bit words (with header+trailer being 16, this remains valid)
0078           if (!size)
0079             size = 16;
0080           randFedSizes[sourceId][eventId] = size;
0081         }
0082       }
0083     }
0084 
0085     for (auto sourceId : sourceIdList_) {
0086       FEDRawData& feddata = data->FEDData(sourceId);
0087 
0088       auto size = sizeof(DTHOrbitHeader_v1);
0089       for (auto eventId : eventIdList_)
0090         size += randFedSizes[sourceId][eventId] + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
0091       feddata.resize(size);
0092 
0093       uint64_t fragments_size_bytes = sizeof(DTHOrbitHeader_v1);
0094       //uint32_t runningChecksum = 0xffffffffU;
0095       uint32_t runningChecksum = 0;
0096       for (auto eventId : eventIdList_) {
0097         unsigned char* fedaddr = feddata.data() + fragments_size_bytes;
0098         //fragments_size_bytes += fillFED(fedaddr, sourceId, eventId, randFedSizes[sourceId][eventId], runningChecksum);
0099         fragments_size_bytes +=
0100             fillSLRFED(fedaddr, sourceId, eventId, orbitId, randFedSizes[sourceId][eventId], runningChecksum);
0101       }
0102       //in place construction
0103       new (feddata.data()) DTHOrbitHeader_v1(sourceId,
0104                                              orbitId,
0105                                              e.id().run(),
0106                                              fragments_size_bytes >> evf::DTH_WORD_NUM_BYTES_SHIFT,
0107                                              eventIdList_.size(),
0108                                              runningChecksum,
0109                                              0);
0110     }
0111   }
0112 
0113   void DTHFakeReader::produce(edm::Event& e, edm::EventSetup const& es) {
0114     edm::Handle<FEDRawDataCollection> rawdata;
0115     FEDRawDataCollection* fedcoll = nullptr;
0116     fillRawData(e, fedcoll);
0117     std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0118     e.put(std::move(bare_product));
0119   }
0120 
0121   uint32_t DTHFakeReader::fillSLRFED(unsigned char* buf,
0122                                      const uint32_t sourceId,
0123                                      edm::EventNumber_t eventId,
0124                                      const uint32_t orbitId,
0125                                      uint32_t size,
0126                                      uint32_t& accum_crc32c) {
0127     // Generate size...
0128     const unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
0129     const unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);
0130 
0131     uint32_t totsize = size + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
0132     const unsigned fragsize = size + h_size_ + t_size_;
0133 
0134     //Fill SLinkRocket header
0135     uint8_t emu_status = 2;  //set 2 indicating fragment generated by DTH
0136     uint16_t l1a_types = 1;  //set provisionally to 1, to be revised later
0137     uint8_t l1a_subtype = 0;
0138     new ((void*)buf) SLinkRocketHeader_v3(eventId, emu_status, l1a_subtype, l1a_types, sourceId);
0139 
0140     // Payload = all 0s or random
0141     if (fillRandom_) {
0142       //fill FED with random values
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*)(buf + h_size_ + i)) = (unsigned int)std::rand();
0146       }
0147       //remainder
0148       for (size_t i = size_ui; i < size; i++) {
0149         *(buf + h_size_ + i) = std::rand() & 0xff;
0150       }
0151     }
0152 
0153     //Fill SLinkRocket trailer
0154     uint16_t crc = 0;  // FIXME : get CRC16
0155     uint16_t bxid = 0;
0156     uint8_t status = 0;
0157     //size is in bytes, it will be converted by constructor
0158     new ((void*)(buf + h_size_ + size))
0159         SLinkRocketTrailer_v3(crc, fragsize >> evf::SLR_WORD_NUM_BYTES_SHIFT, bxid, orbitId, crc, status);
0160 
0161     //fill DTH fragment trailer
0162     void* dthTrailerAddr = buf + fragsize;
0163     new (dthTrailerAddr) DTHFragmentTrailer_v1(fragsize >> evf::DTH_WORD_NUM_BYTES_SHIFT, 0, crc, eventId);
0164 
0165     //accumulate crc32 checksum
0166     accum_crc32c = crc32c(accum_crc32c, (const uint8_t*)buf, totsize);
0167 
0168     return totsize;
0169   }
0170 
0171   uint32_t DTHFakeReader::fillFED(
0172       unsigned char* buf, const int sourceId, edm::EventNumber_t eventId, uint32_t size, uint32_t& accum_crc32c) {
0173     // Generate size...
0174     const unsigned h_size = 8;
0175     const unsigned t_size = 8;
0176 
0177     //header+trailer+payload
0178     uint32_t totsize = size + h_size + t_size + sizeof(DTHFragmentTrailer_v1);
0179 
0180     // Generate header
0181     //FEDHeader::set(feddata.data(),
0182     FEDHeader::set(buf,
0183                    1,          // Trigger type
0184                    eventId,    // LV1_id (24 bits)
0185                    0,          // BX_id
0186                    sourceId);  // source_id
0187 
0188     // Payload = all 0s or random
0189     if (fillRandom_) {
0190       //fill FED with random values
0191       size_t size_ui = size - size % sizeof(unsigned int);
0192       for (size_t i = 0; i < size_ui; i += sizeof(unsigned int)) {
0193         *((unsigned int*)(buf + h_size + i)) = (unsigned int)std::rand();
0194       }
0195       //remainder
0196       for (size_t i = size_ui; i < size; i++) {
0197         *(buf + h_size + i) = std::rand() & 0xff;
0198       }
0199     }
0200 
0201     // Generate trailer
0202     int crc = 0;  // FIXME : get CRC16
0203     FEDTrailer::set(buf + h_size + size,
0204                     size / 8 + 2,  // in 64 bit words
0205                     crc,
0206                     0,   // Evt_stat
0207                     0);  // TTS bits
0208 
0209     //FIXME: accumulate crc32 checksum
0210     //crc32c = 0;
0211 
0212     void* dthTrailerAddr = buf + h_size + t_size + size;
0213     new (dthTrailerAddr)
0214         DTHFragmentTrailer_v1((h_size + t_size + size) >> evf::DTH_WORD_NUM_BYTES_SHIFT, 0, crc, eventId);
0215     return totsize;
0216   }
0217 
0218   void DTHFakeReader::beginLuminosityBlock(edm::LuminosityBlock const& iL, edm::EventSetup const& iE) {
0219     std::cout << "DTHFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
0220     fakeLs_ = iL.luminosityBlock();
0221   }
0222 
0223   void DTHFakeReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0224     edm::ParameterSetDescription desc;
0225     desc.setComment("Injector of generated DTH raw orbit fragments for DSAQ testing");
0226     desc.addUntracked<bool>("fillRandom", false);
0227     desc.addUntracked<unsigned int>("meanSize", 1024);
0228     desc.addUntracked<unsigned int>("width", 1024);
0229     desc.addUntracked<unsigned int>("injectErrPpm", 1024);
0230     desc.addUntracked<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>());
0231     descriptions.add("DTHFakeReader", desc);
0232   }
0233 }  //namespace evf