File indexing completed on 2024-04-06 12:10:53
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <iostream>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/stream/EDFilter.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/Utilities/interface/StreamID.h"
0032
0033 #include "DataFormats/TCDS/interface/TCDSRecord.h"
0034
0035
0036
0037
0038
0039 class TriggerRulePrefireVetoFilter : public edm::stream::EDFilter<> {
0040 public:
0041 explicit TriggerRulePrefireVetoFilter(const edm::ParameterSet&);
0042
0043 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0044
0045 private:
0046 bool filter(edm::Event&, const edm::EventSetup&) override;
0047
0048
0049 edm::EDGetTokenT<TCDSRecord> tcdsRecordToken_;
0050 };
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063 TriggerRulePrefireVetoFilter::TriggerRulePrefireVetoFilter(const edm::ParameterSet& iConfig)
0064 : tcdsRecordToken_(consumes<TCDSRecord>(iConfig.getParameter<edm::InputTag>("tcdsRecordLabel"))) {
0065
0066 }
0067
0068
0069
0070
0071
0072
0073 bool TriggerRulePrefireVetoFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0074 edm::Handle<TCDSRecord> tcdsRecordH;
0075 iEvent.getByToken(tcdsRecordToken_, tcdsRecordH);
0076 const auto& tcdsRecord = *tcdsRecordH.product();
0077
0078 uint64_t thisEvent = (tcdsRecord.getBXID() - 1) + tcdsRecord.getOrbitNr() * 3564ull;
0079
0080 std::vector<uint64_t> eventHistory;
0081 for (auto&& l1a : tcdsRecord.getFullL1aHistory()) {
0082 eventHistory.push_back(thisEvent - ((l1a.getBXID() - 1) + l1a.getOrbitNr() * 3564ull));
0083 }
0084
0085
0086 if (eventHistory.size() < 4) {
0087 edm::LogError("TriggerRulePrefireVetoFilter") << "Unexpectedly small L1A history from TCDSRecord";
0088 }
0089
0090
0091 if (eventHistory[0] < 3ull) {
0092 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (1 in 3)";
0093 }
0094 if (eventHistory[0] == 3ull)
0095 return true;
0096
0097
0098 if (eventHistory[0] < 25ull and eventHistory[1] < 25ull) {
0099 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (2 in 25)";
0100 }
0101 if (eventHistory[0] < 25ull and eventHistory[1] == 25ull)
0102 return true;
0103
0104
0105 if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] < 100ull) {
0106 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (3 in 100)";
0107 }
0108 if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] == 100ull)
0109 return true;
0110
0111
0112 if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] < 240ull) {
0113 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (4 in 240)";
0114 }
0115 if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] == 240ull)
0116 return true;
0117
0118 return false;
0119 }
0120
0121
0122 void TriggerRulePrefireVetoFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0123
0124
0125 edm::ParameterSetDescription desc;
0126 desc.setUnknown();
0127 descriptions.addDefault(desc);
0128 }
0129
0130 DEFINE_FWK_MODULE(TriggerRulePrefireVetoFilter);