Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191
#include "DTHFakeReader.h"
#include "DataFormats/FEDRawData/interface/FEDHeader.h"
#include "DataFormats/FEDRawData/interface/FEDTrailer.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"
#include "DataFormats/FEDRawData/interface/SLinkRocketHeaders.h"
#include "DataFormats/TCDS/interface/TCDSRaw.h"

//#include "EventFilter/Utilities/interface/GlobalEventNumber.h"
#include "EventFilter/Utilities/interface/crc32c.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "CLHEP/Random/RandGauss.h"

#include <cmath>
#include <sys/time.h>
#include <cstring>
#include <cstdlib>
#include <chrono>

//using namespace edm;
namespace evf {

  constexpr unsigned minOrbitBx = 1;
  constexpr unsigned maxOrbitBx = 2464;
  constexpr unsigned avgEventsPerOrbit = 70;

  //constexpr unsigned h_size_ = 8;//for SLink FEDs
  //constexpr unsigned t_size_ = 8;

  constexpr unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
  constexpr unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);

  constexpr double rndFactor = (maxOrbitBx - minOrbitBx + 1) / (double(avgEventsPerOrbit) * RAND_MAX);

  DTHFakeReader::DTHFakeReader(const edm::ParameterSet& pset)
      : fillRandom_(pset.getUntrackedParameter<bool>("fillRandom", false)),
        meansize_(pset.getUntrackedParameter<unsigned int>("meanSize", 1024)),
        width_(pset.getUntrackedParameter<unsigned int>("width", 1024)),
        injected_errors_per_million_events_(pset.getUntrackedParameter<unsigned int>("injectErrPpm", 0)),
        sourceIdList_(
            pset.getUntrackedParameter<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>())),
        modulo_error_events_(injected_errors_per_million_events_ ? 1000000 / injected_errors_per_million_events_
                                                                 : 0xffffffff) {
    if (fillRandom_) {
      //intialize random seed
      auto time_count =
          static_cast<long unsigned int>(std::chrono::high_resolution_clock::now().time_since_epoch().count());
      std::srand(time_count & 0xffffffff);
    }
    produces<RawDataBuffer>();
  }

  void DTHFakeReader::fillRawData(edm::Event& e, RawDataBuffer*& data) {
    // a null pointer is passed, need to allocate the fed collection (reusing it as container)
    //auto ls = e.luminosityBlock();
    //this will be used as orbit counter
    edm::EventNumber_t orbitId = e.id().event();

    //generate eventID. Orbits start from 0 or 1?
    std::vector<uint64_t> eventIdList_;
    std::map<unsigned int, std::map<uint64_t, uint32_t>> randFedSizes;
    for (auto sourceId : sourceIdList_) {
      randFedSizes[sourceId] = std::map<uint64_t, uint32_t>();
    }
    uint32_t totSize = 0;

    //randomize which orbit was accepted
    for (unsigned i = minOrbitBx; i <= maxOrbitBx; i++) {
      if ((std::rand() * rndFactor) < 1) {
        uint64_t eventId = orbitId * maxOrbitBx + i;
        eventIdList_.push_back(eventId);
        for (auto sourceId : sourceIdList_) {
          float logsiz = CLHEP::RandGauss::shoot(std::log(meansize_), std::log(meansize_) - std::log(width_ / 2.));
          size_t size = int(std::exp(logsiz));
          size -= size % 16;  // all blocks aligned to 128 bit words (with header+trailer being 16, this remains valid)
          if (!size)
            size = 16;
          randFedSizes[sourceId][eventId] = size;
        }
      }
    }

    //calculate buffer size and create it
    for (auto sourceId : sourceIdList_) {
      auto size = sizeof(DTHOrbitHeader_v1);
      for (auto eventId : eventIdList_)
        size += randFedSizes[sourceId][eventId] + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
      totSize += size;
    }
    data = new RawDataBuffer(totSize);

    for (auto sourceId : sourceIdList_) {
      auto size = sizeof(DTHOrbitHeader_v1);
      for (auto eventId : eventIdList_)
        size += randFedSizes[sourceId][eventId] + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
      unsigned char* feddata = data->addSource(sourceId, nullptr, size);

      uint64_t fragments_size_bytes = sizeof(DTHOrbitHeader_v1);
      //uint32_t runningChecksum = 0xffffffffU;
      uint32_t runningChecksum = 0;
      for (auto eventId : eventIdList_) {
        unsigned char* fedaddr = feddata + fragments_size_bytes;
        fragments_size_bytes +=
            fillSLRFED(fedaddr, sourceId, eventId, orbitId, randFedSizes[sourceId][eventId], runningChecksum);
      }
      //in place construction
      new (static_cast<void*>(feddata)) DTHOrbitHeader_v1(sourceId,
                                                          e.id().run(),
                                                          orbitId,
                                                          eventIdList_.size(),
                                                          fragments_size_bytes >> evf::DTH_WORD_NUM_BYTES_SHIFT,
                                                          0,
                                                          runningChecksum);
    }
  }

  void DTHFakeReader::produce(edm::Event& e, edm::EventSetup const& es) {
    edm::Handle<RawDataBuffer> rawdata;
    RawDataBuffer* fedcoll = nullptr;
    fillRawData(e, fedcoll);
    std::unique_ptr<RawDataBuffer> bare_product(fedcoll);
    e.put(std::move(bare_product));
  }

  uint32_t DTHFakeReader::fillSLRFED(unsigned char* buf,
                                     const uint32_t sourceId,
                                     edm::EventNumber_t eventId,
                                     const uint32_t orbitId,
                                     uint32_t size,
                                     uint32_t& accum_crc32c) {
    // Generate size...
    const unsigned h_size_ = sizeof(SLinkRocketHeader_v3);
    const unsigned t_size_ = sizeof(SLinkRocketTrailer_v3);

    uint32_t totsize = size + h_size_ + t_size_ + sizeof(DTHFragmentTrailer_v1);
    const unsigned fragsize = size + h_size_ + t_size_;

    //Fill SLinkRocket header
    uint8_t emu_status = 2;  //set 2 indicating fragment generated by DTH (emulator)
    uint16_t l1a_types = 1;  //set provisionally to 1, to be revised later
    uint8_t l1a_phys = 0;
    new ((void*)buf) SLinkRocketHeader_v3(sourceId, l1a_types, l1a_phys, emu_status, eventId);

    // Payload = all 0s or random
    if (fillRandom_) {
      //fill FED with random values
      size_t size_ui = size - size % sizeof(unsigned int);
      for (size_t i = 0; i < size_ui; i += sizeof(unsigned int)) {
        *((unsigned int*)(buf + h_size_ + i)) = (unsigned int)std::rand();
      }
      //remainder
      for (size_t i = size_ui; i < size; i++) {
        *(buf + h_size_ + i) = std::rand() & 0xff;
      }
    }

    //Fill SLinkRocket trailer
    uint16_t crc = 0;  // FIXME : get CRC16
    uint16_t bxid = 0;
    uint8_t status = 0;
    //size is in bytes, it will be converted by constructor
    new ((void*)(buf + h_size_ + size))
        SLinkRocketTrailer_v3(status, crc, orbitId, bxid, fragsize >> evf::SLR_WORD_NUM_BYTES_SHIFT, crc);

    //fill DTH fragment trailer
    void* dthTrailerAddr = buf + fragsize;
    new (dthTrailerAddr) DTHFragmentTrailer_v1(0, fragsize >> evf::DTH_WORD_NUM_BYTES_SHIFT, eventId, crc);

    //accumulate crc32 checksum
    accum_crc32c = crc32c(accum_crc32c, (const uint8_t*)buf, totsize);

    return totsize;
  }

  void DTHFakeReader::beginLuminosityBlock(edm::LuminosityBlock const& iL, edm::EventSetup const& iE) {
    std::cout << "DTHFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
    fakeLs_ = iL.luminosityBlock();
  }

  void DTHFakeReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
    edm::ParameterSetDescription desc;
    desc.setComment("Injector of generated DTH raw orbit fragments for DSAQ testing");
    desc.addUntracked<bool>("fillRandom", false);
    desc.addUntracked<unsigned int>("meanSize", 1024);
    desc.addUntracked<unsigned int>("width", 1024);
    desc.addUntracked<unsigned int>("injectErrPpm", 1024);
    desc.addUntracked<std::vector<unsigned int>>("sourceIdList", std::vector<unsigned int>());
    descriptions.add("DTHFakeReader", desc);
  }
}  //namespace evf