TriggerRulePrefireVetoFilter

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
// -*- C++ -*-
//
// Package:    L1Trigger/PrefireAnalysis
// Class:      TriggerRulePrefireVetoFilter
//
/**\class TriggerRulePrefireVetoFilter TriggerRulePrefireVetoFilter.cc L1Trigger/PrefireAnalysis/plugins/TriggerRulePrefireVetoFilter.cc

 Description: [one line class summary]

 Implementation:
    [Notes on implementation]
*/
//
// Original Author:  Nicholas Charles Smith
//         Created:  Mon, 28 May 2018 15:39:05 GMT
//
//

// system include files
#include <memory>
#include <iostream>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/stream/EDFilter.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"

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

#include "DataFormats/TCDS/interface/TCDSRecord.h"

//
// class declaration
//

class TriggerRulePrefireVetoFilter : public edm::stream::EDFilter<> {
public:
  explicit TriggerRulePrefireVetoFilter(const edm::ParameterSet&);

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
  bool filter(edm::Event&, const edm::EventSetup&) override;

  // ----------member data ---------------------------
  edm::EDGetTokenT<TCDSRecord> tcdsRecordToken_;
};

//
// constants, enums and typedefs
//

//
// static data member definitions
//

//
// constructors and destructor
//
TriggerRulePrefireVetoFilter::TriggerRulePrefireVetoFilter(const edm::ParameterSet& iConfig)
    : tcdsRecordToken_(consumes<TCDSRecord>(iConfig.getParameter<edm::InputTag>("tcdsRecordLabel"))) {
  //now do what ever initialization is needed
}

//
// member functions
//

// ------------ method called on each new Event  ------------
bool TriggerRulePrefireVetoFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
  edm::Handle<TCDSRecord> tcdsRecordH;
  iEvent.getByToken(tcdsRecordToken_, tcdsRecordH);
  const auto& tcdsRecord = *tcdsRecordH.product();

  uint64_t thisEvent = (tcdsRecord.getBXID() - 1) + tcdsRecord.getOrbitNr() * 3564ull;

  std::vector<uint64_t> eventHistory;
  for (auto&& l1a : tcdsRecord.getFullL1aHistory()) {
    eventHistory.push_back(thisEvent - ((l1a.getBXID() - 1) + l1a.getOrbitNr() * 3564ull));
  }

  // should be 16 according to TCDSRecord.h, we only care about the last 4
  if (eventHistory.size() < 4) {
    edm::LogError("TriggerRulePrefireVetoFilter") << "Unexpectedly small L1A history from TCDSRecord";
  }

  // No more than 1 L1A in 3 BX
  if (eventHistory[0] < 3ull) {
    edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (1 in 3)";
  }
  if (eventHistory[0] == 3ull)
    return true;

  // No more than 2 L1As in 25 BX
  if (eventHistory[0] < 25ull and eventHistory[1] < 25ull) {
    edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (2 in 25)";
  }
  if (eventHistory[0] < 25ull and eventHistory[1] == 25ull)
    return true;

  // No more than 3 L1As in 100 BX
  if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] < 100ull) {
    edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (3 in 100)";
  }
  if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] == 100ull)
    return true;

  // No more than 4 L1As in 240 BX
  if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] < 240ull) {
    edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (4 in 240)";
  }
  if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] == 240ull)
    return true;

  return false;
}

// ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
void TriggerRulePrefireVetoFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  //The following says we do not know what parameters are allowed so do no validation
  // Please change this to state exactly what you do use, even if it is no parameters
  edm::ParameterSetDescription desc;
  desc.setUnknown();
  descriptions.addDefault(desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(TriggerRulePrefireVetoFilter);