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
#include "LaserAlignmentEventFilter.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/EventSetupRecordKey.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "CondFormats/SiStripObjects/interface/SiStripFedCabling.h"
#include "EventFilter/SiStripRawToDigi/interface/SiStripFEDBuffer.h"

#include <algorithm>

// make a sorted copy of a vector, optionally of a different type
template <typename T, typename U>
std::vector<T> copy_and_sort_vector(std::vector<U> const& input) {
  std::vector<T> copy(input.begin(), input.end());
  std::sort(copy.begin(), copy.end());
  return copy;
}

///
/// constructors and destructor
///
LaserAlignmentEventFilter::LaserAlignmentEventFilter(const edm::ParameterSet& iConfig)
    : cablingToken_(esConsumes()),
      FED_collection_token(consumes<FEDRawDataCollection>(iConfig.getParameter<edm::InputTag>("FedInputTag"))),
      las_fed_ids(copy_and_sort_vector<uint16_t>(iConfig.getParameter<std::vector<int>>("FED_IDs"))),
      las_signal_ids(copy_and_sort_vector<uint32_t>(iConfig.getParameter<std::vector<int>>("SIGNAL_IDs"))),
      single_channel_thresh(iConfig.getParameter<unsigned>("SINGLE_CHANNEL_THRESH")),
      channel_count_thresh(iConfig.getParameter<unsigned>("CHANNEL_COUNT_THRESH")) {}

///
///
///
LaserAlignmentEventFilter::~LaserAlignmentEventFilter() {}

/// Checks for Laser Signals in specific modules
/// For Accessing FED Data see also EventFilter/SiStripRawToDigi/src/SiStripRawToDigiUnpacker.cc and related files
bool LaserAlignmentEventFilter::filter(edm::StreamID sid, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
  //unsigned int  det_ctr = 0;      // count how many modules are tested for signal
  unsigned int sig_ctr = 0;  // count how many modules have signal
  //unsigned long buffer_sum = 0;   // sum of buffer sizes

  // retrieve FED raw data (by label, which is "source" by default)
  edm::Handle<FEDRawDataCollection> buffers;
  iEvent.getByToken(FED_collection_token, buffers);

  // read the cabling map from the EventSetup
  const auto& cabling = &iSetup.getData(cablingToken_);

  std::vector<uint16_t>::const_iterator ifed = las_fed_ids.begin();
  for (; ifed != las_fed_ids.end(); ifed++) {
    // Retrieve FED raw data for given FED
    const FEDRawData& input = buffers->FEDData(static_cast<int>(*ifed));
    LogDebug("LaserAlignmentEventFilter") << "Examining FED " << *ifed;

    // check on FEDRawData pointer
    if (not input.data())
      continue;

    // check on FEDRawData size
    if (not input.size())
      continue;

    // construct FEDBuffer
    const auto st_buffer = sistrip::preconstructCheckFEDBuffer(input);
    if (sistrip::FEDBufferStatusCode::SUCCESS != st_buffer) {
      throw cms::Exception("FEDBuffer") << st_buffer << " (check debug output for more details)";
    }
    sistrip::FEDBuffer buffer{input};
    const auto st_chan = buffer.findChannels();
    if (sistrip::FEDBufferStatusCode::SUCCESS != st_chan) {
      throw cms::Exception("FEDBuffer") << st_chan << " (check debug output for more details)";
    }
    if (not buffer.doChecks(true)) {
      edm::LogWarning("LaserAlignmentEventFilter") << "FED Buffer check fails for FED ID " << *ifed << ".";
      continue;
    }

    // get the cabling connections for this FED
    auto const& conns = cabling->fedConnections(*ifed);

    // Iterate through FED channels, extract payload and create Digis
    std::vector<FedChannelConnection>::const_iterator iconn = conns.begin();
    for (; iconn != conns.end(); iconn++) {
      if (std::binary_search(las_signal_ids.begin(), las_signal_ids.end(), iconn->detId())) {
        LogDebug("LaserAlignmentEventFilter")
            << " Found LAS signal module in FED " << *ifed << "  DetId: " << iconn->detId() << "\n"
            << "buffer.channel(iconn->fedCh()).size(): " << buffer.channel(iconn->fedCh()).length();
        //++det_ctr;
        if (buffer.channel(iconn->fedCh()).length() > single_channel_thresh)
          ++sig_ctr;
        //buffer_sum += buffer.channel(iconn->fedCh()).length();
        if (sig_ctr > channel_count_thresh) {
          LogDebug("LaserAlignmentEventFilter") << "Event identified as LAS";
          return true;
        }
      }
    }  // channel loop
  }  // FED loop

  //   LogDebug("LaserAlignmentEventFilter") << det_ctr << " channels were tested for signal\n"
  // 			    <<sig_ctr << " channels have signal\n"
  // 			    << "Sum of buffer sizes: " << buffer_sum;

  return false;
}

//define this as a plug-in
DEFINE_FWK_MODULE(LaserAlignmentEventFilter);