Skip to content

Commit 8be9981

Browse files
author
Benjamin Huber
committed
Add L1GTAcceptFilter
1 parent 6fc759d commit 8be9981

File tree

2 files changed

+68
-0
lines changed

2 files changed

+68
-0
lines changed
Lines changed: 62 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,62 @@
1+
#include "FWCore/Framework/interface/Event.h"
2+
#include "FWCore/Framework/interface/global/EDFilter.h"
3+
#include "FWCore/ParameterSet/interface/ParameterSet.h"
4+
#include "DataFormats/L1Trigger/interface/P2GTAlgoBlock.h"
5+
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
6+
7+
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
8+
#include "FWCore/ParameterSet/interface/allowedValues.h"
9+
#include "FWCore/Framework/interface/MakerMacros.h"
10+
#include "DataFormats/Common/interface/Handle.h"
11+
#include "FWCore/Utilities/interface/EDGetToken.h"
12+
13+
using namespace l1t;
14+
15+
class L1GTAcceptFilter : public edm::global::EDFilter<> {
16+
public:
17+
explicit L1GTAcceptFilter(const edm::ParameterSet&);
18+
~L1GTAcceptFilter() override = default;
19+
20+
static void fillDescriptions(edm::ConfigurationDescriptions&);
21+
22+
enum DecisionType { beforeBxMaskAndPrescale, beforePrescale, final };
23+
24+
private:
25+
bool filter(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
26+
27+
const edm::EDGetTokenT<P2GTAlgoBlockMap> algoBlocksToken_;
28+
const DecisionType decisionEnum_;
29+
};
30+
31+
L1GTAcceptFilter::L1GTAcceptFilter(const edm::ParameterSet& config)
32+
: algoBlocksToken_(consumes<P2GTAlgoBlockMap>(config.getParameter<edm::InputTag>("algoBlocksTag"))),
33+
decisionEnum_(config.getParameter<std::string>("decision") == "beforeBxMaskAndPrescale" ? beforeBxMaskAndPrescale
34+
: config.getParameter<std::string>("decision") == "beforePrescale" ? beforePrescale
35+
: final) {}
36+
37+
bool L1GTAcceptFilter::filter(edm::StreamID, edm::Event& event, const edm::EventSetup& setup) const {
38+
const P2GTAlgoBlockMap& algoMap = event.get(algoBlocksToken_);
39+
bool decision = false;
40+
for (const auto& [name, algoBlock] : algoMap) {
41+
if (decisionEnum_ == beforeBxMaskAndPrescale) {
42+
decision |= algoBlock.decisionBeforeBxMaskAndPrescale();
43+
} else if (decisionEnum_ == beforePrescale) {
44+
decision |= algoBlock.decisionBeforePrescale();
45+
} else {
46+
decision |= algoBlock.decisionFinal();
47+
}
48+
}
49+
50+
return decision;
51+
}
52+
53+
void L1GTAcceptFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
54+
edm::ParameterSetDescription desc;
55+
desc.add<edm::InputTag>("algoBlocksTag");
56+
desc.ifValue(edm::ParameterDescription<std::string>("decision", "final", true),
57+
edm::allowedValues<std::string>("beforeBxMaskAndPrescale", "beforePrescale", "final"));
58+
59+
descriptions.addWithDefaultLabel(desc);
60+
}
61+
62+
DEFINE_FWK_MODULE(L1GTAcceptFilter);
Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,6 @@
1+
import FWCore.ParameterSet.Config as cms
2+
3+
l1tGTAcceptFilter = cms.EDFilter(
4+
"L1GTAcceptFilter",
5+
algoBlocksTag = cms.InputTag("l1tGTAlgoBlockProducer"),
6+
)

0 commit comments

Comments
 (0)