File indexing completed on 2023-10-25 09:45:39
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 ~TriggerRulePrefireVetoFilter() override;
0043
0044 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0045
0046 private:
0047 void beginStream(edm::StreamID) override;
0048 bool filter(edm::Event&, const edm::EventSetup&) override;
0049 void endStream() override;
0050
0051
0052
0053
0054
0055
0056
0057 edm::EDGetTokenT<TCDSRecord> tcdsRecordToken_;
0058 };
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071 TriggerRulePrefireVetoFilter::TriggerRulePrefireVetoFilter(const edm::ParameterSet& iConfig)
0072 : tcdsRecordToken_(consumes<TCDSRecord>(iConfig.getParameter<edm::InputTag>("tcdsRecordLabel"))) {
0073
0074 }
0075
0076 TriggerRulePrefireVetoFilter::~TriggerRulePrefireVetoFilter() {
0077
0078
0079 }
0080
0081
0082
0083
0084
0085
0086 bool TriggerRulePrefireVetoFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0087 edm::Handle<TCDSRecord> tcdsRecordH;
0088 iEvent.getByToken(tcdsRecordToken_, tcdsRecordH);
0089 const auto& tcdsRecord = *tcdsRecordH.product();
0090
0091 uint64_t thisEvent = (tcdsRecord.getBXID() - 1) + tcdsRecord.getOrbitNr() * 3564ull;
0092
0093 std::vector<uint64_t> eventHistory;
0094 for (auto&& l1a : tcdsRecord.getFullL1aHistory()) {
0095 eventHistory.push_back(thisEvent - ((l1a.getBXID() - 1) + l1a.getOrbitNr() * 3564ull));
0096 }
0097
0098
0099 if (eventHistory.size() < 4) {
0100 edm::LogError("TriggerRulePrefireVetoFilter") << "Unexpectedly small L1A history from TCDSRecord";
0101 }
0102
0103
0104 if (eventHistory[0] < 3ull) {
0105 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (1 in 3)";
0106 }
0107 if (eventHistory[0] == 3ull)
0108 return true;
0109
0110
0111 if (eventHistory[0] < 25ull and eventHistory[1] < 25ull) {
0112 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (2 in 25)";
0113 }
0114 if (eventHistory[0] < 25ull and eventHistory[1] == 25ull)
0115 return true;
0116
0117
0118 if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] < 100ull) {
0119 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (3 in 100)";
0120 }
0121 if (eventHistory[0] < 100ull and eventHistory[1] < 100ull and eventHistory[2] == 100ull)
0122 return true;
0123
0124
0125 if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] < 240ull) {
0126 edm::LogError("TriggerRulePrefireVetoFilter") << "Found an L1A in an impossible location?! (4 in 240)";
0127 }
0128 if (eventHistory[0] < 240ull and eventHistory[1] < 240ull and eventHistory[2] < 240ull and eventHistory[3] == 240ull)
0129 return true;
0130
0131 return false;
0132 }
0133
0134
0135 void TriggerRulePrefireVetoFilter::beginStream(edm::StreamID) {}
0136
0137
0138 void TriggerRulePrefireVetoFilter::endStream() {}
0139
0140
0141
0142
0143
0144
0145
0146
0147
0148
0149
0150
0151
0152
0153
0154
0155
0156
0157
0158
0159
0160
0161
0162
0163
0164
0165
0166
0167
0168
0169
0170
0171
0172
0173 void TriggerRulePrefireVetoFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0174
0175
0176 edm::ParameterSetDescription desc;
0177 desc.setUnknown();
0178 descriptions.addDefault(desc);
0179 }
0180
0181 DEFINE_FWK_MODULE(TriggerRulePrefireVetoFilter);